Octopus
geom_opt.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2007 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module geom_opt_oct_m
22 use accel_oct_m
24 use debug_oct_m
28 use global_oct_m
30 use io_oct_m
33 use ions_oct_m
34 use, intrinsic :: iso_fortran_env
36 use lcao_oct_m
37 use loct_oct_m
38 use math_oct_m
39 use mesh_oct_m
42 use mpi_oct_m
45 use parser_oct_m
46 use pcm_oct_m
50 use scf_oct_m
56 use unit_oct_m
58 use v_ks_oct_m
60
61 implicit none
62
63 private
64 public :: geom_opt_run
65
66 type geom_opt_t
67 private
68 integer(int64) :: type
69 integer :: method
70 real(real64) :: step
71 real(real64) :: line_tol
72 real(real64) :: fire_mass
73 integer :: fire_integrator
74 real(real64) :: tolgrad
75 real(real64) :: toldr
76 integer :: max_iter
77 integer :: what2minimize
78
80 type(scf_t) :: scfv
81 type(ions_t), pointer :: ions
82 type(hamiltonian_elec_t), pointer :: hm
83 type(electrons_t), pointer :: syst
84 class(mesh_t), pointer :: mesh
85 type(states_elec_t), pointer :: st
86 integer :: dim
87 integer :: periodic_dim
88 integer :: size
89 integer :: fixed_atom = 0
90
91 real(real64), allocatable :: cell_force(:, :)
92 logical :: symmetrize = .false.
93 real(real64), allocatable :: initial_length(:)
94 real(real64), allocatable :: initial_rlattice(:, :)
95 real(real64), allocatable :: inv_initial_rlattice(:, :)
96 real(real64) :: pressure
97 end type geom_opt_t
98
99 type(geom_opt_t), save :: g_opt
100
101 integer, parameter :: &
102 MINWHAT_ENERGY = 1, &
104
105 integer, parameter :: &
106 GO_IONS = 1, &
107 go_cell = 2, &
108 go_volume = 4
109
110contains
111
112 ! ---------------------------------------------------------
113 subroutine geom_opt_run(system, from_scratch)
114 class(*), intent(inout) :: system
115 logical, intent(inout) :: from_scratch
116
117 push_sub(geom_opt_run)
118
119 select type (system)
120 class is (multisystem_basic_t)
121 message(1) = "CalculationMode = go not implemented for multi-system calculations"
122 call messages_fatal(1, namespace=system%namespace)
123 type is (electrons_t)
124 call geom_opt_run_legacy(system, from_scratch)
125 end select
126
127 pop_sub(geom_opt_run)
128 end subroutine geom_opt_run
129
130 ! ---------------------------------------------------------
131 subroutine geom_opt_run_legacy(sys, fromscratch)
132 type(electrons_t), target, intent(inout) :: sys
133 logical, intent(inout) :: fromscratch
134
135 integer :: ierr
136 real(real64), allocatable :: coords(:)
137 real(real64) :: energy
138
139 real(real64), allocatable :: mass(:)
140 integer :: iatom, imass
141 type(restart_t) :: restart_load
142
143 push_sub(geom_opt_run_legacy)
144
145 if (sys%space%periodic_dim == 1 .or. sys%space%periodic_dim == 2) then
146 message(1) = "Geometry optimization not allowed for systems periodic in 1D and 2D, "
147 message(2) = "as in those cases the total energy is not correct."
148 call messages_fatal(2, namespace=sys%namespace)
149 end if
150
151
152 if (sys%hm%pcm%run_pcm) then
153 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
154 end if
155
156 if (sys%kpoints%use_symmetries) then
157 call messages_experimental("KPoints symmetries with CalculationMode = go", namespace=sys%namespace)
158 end if
160 call init_(fromscratch)
162 ! load wavefunctions
163 if (.not. fromscratch) then
164 call restart_init(restart_load, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr)
165 if (ierr == 0) then
166 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
167 end if
168 call restart_end(restart_load)
169 if (ierr /= 0) then
170 message(1) = "Unable to read wavefunctions: Starting from scratch."
171 call messages_warning(1, namespace=sys%namespace)
172 fromscratch = .true.
173 end if
174 end if
176 call scf_init(g_opt%scfv, sys%namespace, sys%gr, sys%ions, sys%st, sys%mc, sys%hm, sys%space)
178 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0) then
179 if (.not. g_opt%scfv%calc_stress) then
180 message(1) = "In order to optimize the cell, one needs to set SCFCalculateStress = yes."
181 call messages_fatal(1, namespace=sys%namespace)
182 end if
183 end if
185 if (fromscratch) then
186 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm, lmm_r = g_opt%scfv%lmm_r)
187 else
188 ! setup Hamiltonian
189 message(1) = 'Info: Setting up Hamiltonian.'
190 call messages_info(1, namespace=sys%namespace)
191 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
192 end if
193
194 g_opt%symmetrize = sys%kpoints%use_symmetries .or. sys%st%symmetrize_density
195
196 !Initial point
197 safe_allocate(coords(1:g_opt%size))
198 call to_coords(g_opt, coords)
199
200 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
201
202 !Minimize
203 select case (g_opt%method)
205 call minimize_multidim_nograd(g_opt%method, g_opt%size, coords, real(g_opt%step, real64) ,&
206 real(g_opt%toldr, real64) , g_opt%max_iter, &
207 calc_point_ng, write_iter_info_ng, energy, ierr)
208 case (minmethod_fire)
209
210 safe_allocate(mass(1:g_opt%size))
211 mass = g_opt%fire_mass
212 imass = 1
213 do iatom = 1, sys%ions%natoms
214 if (g_opt%fixed_atom == iatom) cycle
215 if (g_opt%ions%fixed(iatom)) cycle
216 if (g_opt%fire_mass <= m_zero) mass(imass:imass + 2) = sys%ions%mass(iatom)
217 imass = imass + 3
218 end do
219
220 !TODO: add variable to use Euler integrator
221 call minimize_fire(g_opt%size, g_opt%ions%space%dim, coords, real(g_opt%step, real64) , real(g_opt%tolgrad, real64) , &
222 g_opt%max_iter, calc_point, write_iter_info, energy, ierr, mass, integrator=g_opt%fire_integrator)
223 safe_deallocate_a(mass)
225 case default
226 call minimize_multidim(g_opt%method, g_opt%size, coords, real(g_opt%step, real64) ,&
227 real(g_opt%line_tol, real64) , real(g_opt%tolgrad, real64) , real(g_opt%toldr, real64) , g_opt%max_iter, &
228 calc_point, write_iter_info, energy, ierr)
229 end select
230
231 if (ierr == 1025) then
232 ! not a GSL error, set by our minimize routines, so we must handle it separately
233 message(1) = "Reached maximum number of iterations allowed by GOMaxIter."
234 call messages_info(1, namespace=sys%namespace)
235 else if (ierr /= 0 .and. g_opt%method /= minmethod_fire) then
236 message(1) = "Error occurred during the GSL minimization procedure:"
237 call loct_strerror(ierr, message(2))
238 call messages_fatal(2, namespace=sys%namespace)
239 end if
240
241 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%unpack()
242
243
244 ! print out geometry
245 message(1) = "Writing final coordinates to min.xyz"
246 call messages_info(1, namespace=sys%namespace)
247 call from_coords(g_opt, coords)
248 call g_opt%ions%write_xyz('./min')
249
250 safe_deallocate_a(coords)
251 call scf_end(g_opt%scfv)
252 ! Because g_opt has the "save" atribute, we need to explicitly empty the criteria list here, or there will be a memory leak.
253 call g_opt%scfv%criterion_list%empty()
254 call end_()
255
256 pop_sub(geom_opt_run_legacy)
257 contains
258
259 ! ---------------------------------------------------------
260 subroutine init_(fromscratch)
261 logical, intent(inout) :: fromscratch
262
263 logical :: center, does_exist
264 integer :: iter, iatom, idir
265 character(len=100) :: filename
266 real(real64) :: default_toldr
267 real(real64) :: default_step
268 type(read_coords_info) :: xyz
269
270 push_sub(geom_opt_run_legacy.init_)
271
272 if (sys%space%is_periodic()) then
273 call messages_experimental('Geometry optimization for periodic systems', namespace=sys%namespace)
274 end if
275
276 !%Variable GOType
277 !%Type flag
278 !%Default ions
279 !%Section Calculation Modes::Geometry Optimization
280 !%Description
281 !% This variable defines which parameters are allowed to change during the optimization.
282 !% Multiple options can be chosen e.g. “ion_positions + cell_shape”.
283 !% Only one type of lattice vectors relaxation is possible.
284 !%Option ion_positions 1
285 !% Relax position of ions based on the forces acting on the ions.
286 !%Option cell_shape 2
287 !% Relax cell shape. This changes lattice vector lengths and directions
288 !% based on the stress acting on the lattice vectors.
289 !% See for instance Wentzcovitch, PRB 44, 2358 (1991).
290 !%Option cell_volume 4
291 !% Relax cell volume. Only allow for rescaling the lengths of lattice vectors.
292 !% This is a simplication of the option cell_shape, where only a diagonal strain is allowed.
293 !%End
294
295 call parse_variable(sys%namespace, 'GOType', go_ions, g_opt%type)
296 if (.not. varinfo_valid_option('GOType', g_opt%type, is_flag=.true.)) then
297 call messages_input_error(sys%namespace, 'GOType')
298 end if
299
300 write(message(1),'(a)') 'Input: [GOType = '
301 if (bitand(g_opt%type, go_ions) /= 0) then
302 write(message(1),'(a,1x,a)') trim(message(1)), 'ion_positions'
303 end if
304 if (bitand(g_opt%type, go_cell) /= 0) then
305 if (len_trim(message(1)) > 16) then
306 write(message(1),'(a,1x,a)') trim(message(1)), '+'
307 end if
308 write(message(1),'(a,1x,a)') trim(message(1)), 'cell_shape'
309 end if
310 if (bitand(g_opt%type, go_volume) /= 0) then
311 if (len_trim(message(1)) > 16) then
312 write(message(1),'(a,1x,a)') trim(message(1)), '+'
313 end if
314 write(message(1),'(a,1x,a)') trim(message(1)), 'cell_volume'
315 end if
316 write(message(1),'(2a)') trim(message(1)), ']'
317 call messages_info(1, namespace=sys%namespace)
318
319 if (bitand(g_opt%type, go_volume) /= 0 .and. bitand(g_opt%type, go_cell) /= 0) then
320 message(1) = "Cell and volume optimization cannot be used simultaneously."
321 call messages_fatal(1, namespace=sys%namespace)
322 end if
323
324
325 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0) then
326 if (parse_is_defined(sys%namespace, 'TDMomentumTransfer') .or. &
327 parse_is_defined(sys%namespace, 'TDReducedMomentumTransfer')) then
328 call messages_not_implemented("Cell dynamics with TDMomentumTransfer and TDReducedMomentumTransfer")
329 end if
330 end if
331
332 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0) then
333 if (accel_is_enabled()) then
334 message(1) = "Cell dynamics not supported on GPUs."
335 call messages_fatal(1, namespace=sys%namespace)
336 end if
337 end if
338
339
340 do iatom = 1, sys%ions%natoms
341 select type(spec=>sys%ions%atom(iatom)%species)
342 class is(allelectron_t)
343 write(message(1),'(a)') "Geometry optimization for all-electron potential is not implemented."
344 call messages_fatal(1)
345 end select
346 end do
347
348
349 call states_elec_allocate_wfns(sys%st, sys%gr, packed=.true.)
350
351 ! shortcuts
352 g_opt%mesh => sys%gr
353 g_opt%ions => sys%ions
354 g_opt%st => sys%st
355 g_opt%hm => sys%hm
356 g_opt%syst => sys
357 g_opt%dim = sys%space%dim
358 g_opt%periodic_dim = sys%space%periodic_dim
359
360 g_opt%size = 0
361 ! Ion dyamics
362 if (bitand(g_opt%type, go_ions) /= 0) then
363 g_opt%size = g_opt%dim * g_opt%ions%natoms
364 end if
365
366 ! Cell dynamics
367 if (bitand(g_opt%type, go_cell) /= 0) then
368 g_opt%size = g_opt%size + (g_opt%periodic_dim +1) * g_opt%periodic_dim / 2
369 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
370 end if
371
372 ! Volume dynamics
373 if (bitand(g_opt%type, go_volume) /= 0) then
374 g_opt%size = g_opt%size + g_opt%periodic_dim
375 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:1))
376 ! Store the length of the original lattic vectors, to work with reduced lengthes
377 safe_allocate(g_opt%initial_length(1:g_opt%periodic_dim))
378 do idir = 1, g_opt%periodic_dim
379 g_opt%initial_length(idir) = norm2(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))
380 end do
381 end if
382
383 ! Store the initial lattice vectors and the inverse matrix
384 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0) then
385 safe_allocate(g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
386 g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim) &
387 = g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
388 safe_allocate(g_opt%inv_initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
389 g_opt%inv_initial_rlattice(:, :) = g_opt%initial_rlattice(:, :)
390 call lalg_inverse(g_opt%periodic_dim, g_opt%inv_initial_rlattice, 'dir')
391 end if
392
393 if(g_opt%ions%space%is_periodic()) then
394 call parse_variable(sys%namespace, 'HydrostaticPressure', m_zero, g_opt%pressure)
395 end if
396
397 assert(g_opt%size > 0)
398
399 !%Variable GOCenter
400 !%Type logical
401 !%Default no
402 !%Section Calculation Modes::Geometry Optimization
403 !%Description
404 !% (Experimental) If set to yes, Octopus centers the geometry at
405 !% every optimization step. It also reduces the degrees of
406 !% freedom of the optimization by using the translational
407 !% invariance.
408 !%End
409 call parse_variable(sys%namespace, 'GOCenter', .false., center)
410
411 if (center) then
412 g_opt%fixed_atom = 1
413 g_opt%size = g_opt%size - g_opt%dim
414 call messages_experimental('GOCenter', namespace=sys%namespace)
415 end if
416
417 !Check if atoms are allowed to move and redifine g_opt%size
418 do iatom = 1, g_opt%ions%natoms
419 if (g_opt%ions%fixed(iatom)) then
420 g_opt%size = g_opt%size - g_opt%dim
421 end if
422 end do
423
424 !%Variable GOMethod
425 !%Type integer
426 !%Default fire
427 !%Section Calculation Modes::Geometry Optimization
428 !%Description
429 !% Method by which the minimization is performed. For more information see the
430 !% <a href=http://www.gnu.org/software/gsl/manual/html_node/Multidimensional-Minimization.html>
431 !% GSL documentation</a>.
432 !%Option steep 1
433 !% Simple steepest descent.
434 !%Option steep_native -1
435 !% (Experimental) Non-gsl implementation of steepest descent.
436 !%Option cg_fr 2
437 !% Fletcher-Reeves conjugate-gradient algorithm. The
438 !% conjugate-gradient algorithm proceeds as a succession of line
439 !% minimizations. The sequence of search directions is used to build
440 !% up an approximation to the curvature of the function in the
441 !% neighborhood of the minimum.
442 !%Option cg_pr 3
443 !% Polak-Ribiere conjugate-gradient algorithm.
444 !%Option cg_bfgs 4
445 !% Vector Broyden-Fletcher-Goldfarb-Shanno (BFGS) conjugate-gradient algorithm.
446 !% It is a quasi-Newton method which builds up an approximation to the second
447 !% derivatives of the function <i>f</i> using the difference between successive gradient
448 !% vectors. By combining the first and second derivatives, the algorithm is able
449 !% to take Newton-type steps towards the function minimum, assuming quadratic
450 !% behavior in that region.
451 !%Option cg_bfgs2 5
452 !% The bfgs2 version of this minimizer is the most efficient version available,
453 !% and is a faithful implementation of the line minimization scheme described in
454 !% Fletcher, <i>Practical Methods of Optimization</i>, Algorithms 2.6.2 and 2.6.4.
455 !%Option simplex 6
456 !% This is experimental, and in fact, <b>not</b> recommended unless you just want to
457 !% fool around. It is the Nead-Melder simplex algorithm, as implemented in the
458 !% GNU Scientific Library (GSL). It does not make use of the gradients (<i>i.e.</i>, the
459 !% forces) which makes it less efficient than other schemes. It is included here
460 !% for completeness, since it is free.
461 !%Option fire 8
462 !% The FIRE algorithm. See also <tt>GOFireMass</tt> and <tt>GOFireIntegrator</tt>.
463 !% Ref: E. Bitzek, P. Koskinen, F. Gahler, M. Moseler, and P. Gumbsch, <i>Phys. Rev. Lett.</i> <b>97</b>, 170201 (2006).
464 !%End
465 call parse_variable(sys%namespace, 'GOMethod', minmethod_fire, g_opt%method)
466 if (.not. varinfo_valid_option('GOMethod', g_opt%method)) call messages_input_error(sys%namespace, 'GOMethod')
467
468
469 call messages_print_var_option("GOMethod", g_opt%method, namespace=sys%namespace)
470
471 !%Variable GOTolerance
472 !%Type float
473 !%Default 0.001 H/b (0.051 eV/A)
474 !%Section Calculation Modes::Geometry Optimization
475 !%Description
476 !% Convergence criterion, for stopping the minimization. In
477 !% units of force; minimization is stopped when all forces on
478 !% ions are smaller than this criterion, or the
479 !% <tt>GOMinimumMove</tt> is satisfied. If <tt>GOTolerance < 0</tt>,
480 !% this criterion is ignored.
481 !%End
482 call parse_variable(sys%namespace, 'GOTolerance', 0.001_real64, g_opt%tolgrad, units_inp%force)
483
484 !%Variable GOMinimumMove
485 !%Type float
486 !%Section Calculation Modes::Geometry Optimization
487 !%Description
488 !% Convergence criterion, for stopping the minimization. In
489 !% units of length; minimization is stopped when the coordinates
490 !% of all species change less than <tt>GOMinimumMove</tt>, or the
491 !% <tt>GOTolerance</tt> criterion is satisfied.
492 !% If <tt>GOMinimumMove < 0</tt>, this criterion is ignored.
493 !% Default is -1, except 0.001 b with <tt>GOMethod = simplex</tt>.
494 !% Note that if you use <tt>GOMethod = simplex</tt>,
495 !% then you must supply a non-zero <tt>GOMinimumMove</tt>.
496 !%End
497 if (g_opt%method == minmethod_nmsimplex) then
498 default_toldr = 0.001_real64
499 else
500 default_toldr = -m_one
501 end if
502 call parse_variable(sys%namespace, 'GOMinimumMove', default_toldr, g_opt%toldr, units_inp%length)
503
504 if (g_opt%method == minmethod_nmsimplex .and. g_opt%toldr <= m_zero) call messages_input_error(sys%namespace, 'GOMinimumMove')
505
506 !%Variable GOStep
507 !%Type float
508 !%Section Calculation Modes::Geometry Optimization
509 !%Description
510 !% Initial step for the geometry optimizer. The default is 0.5.
511 !% WARNING: in some weird units.
512 !% For the FIRE minimizer, default value is 0.1 fs,
513 !% and corresponds to the initial time-step for the MD.
514 !%End
515 if (g_opt%method /= minmethod_fire) then
516 default_step = m_half
517 call parse_variable(sys%namespace, 'GOStep', default_step, g_opt%step)
518 else
519 default_step = 0.1_real64*unit_femtosecond%factor
520 call parse_variable(sys%namespace, 'GOStep', default_step, g_opt%step, unit = units_inp%time)
521 end if
522
523 !%Variable GOLineTol
524 !%Type float
525 !%Default 0.1
526 !%Section Calculation Modes::Geometry Optimization
527 !%Description
528 !% Tolerance for line-minimization. Applies only to GSL methods
529 !% that use the forces.
530 !% WARNING: in some weird units.
531 !%End
532 call parse_variable(sys%namespace, 'GOLineTol', 0.1_real64, g_opt%line_tol)
533
534 !%Variable GOMaxIter
535 !%Type integer
536 !%Default 200
537 !%Section Calculation Modes::Geometry Optimization
538 !%Description
539 !% Even if the convergence criterion is not satisfied, the minimization will stop
540 !% after this number of iterations.
541 !%End
542 call parse_variable(sys%namespace, 'GOMaxIter', 200, g_opt%max_iter)
543 if (g_opt%max_iter <= 0) then
544 message(1) = "GOMaxIter has to be larger than 0"
545 call messages_fatal(1, namespace=sys%namespace)
546 end if
547
548 !%Variable GOFireMass
549 !%Type float
550 !%Default 1.0 amu
551 !%Section Calculation Modes::Geometry Optimization
552 !%Description
553 !% The Fire algorithm (<tt>GOMethod = fire</tt>) assumes that all degrees of freedom
554 !% are comparable. All the velocities should be on the same
555 !% scale, which for heteronuclear systems can be roughly
556 !% achieved by setting all the atom masses equal, to the value
557 !% specified by this variable.
558 !% By default the mass of a proton is selected (1 amu).
559 !% However, a selection of <tt>GOFireMass = 0.01</tt> can, in manys systems,
560 !% speed up the geometry optimization procedure.
561 !% If <tt>GOFireMass</tt> <= 0, the masses of each
562 !% species will be used.
563 !%End
564 call parse_variable(sys%namespace, 'GOFireMass', m_one*unit_amu%factor, g_opt%fire_mass, unit = unit_amu)
565
566 !%Variable GOFireIntegrator
567 !%Type integer
568 !%Default verlet
569 !%Section Calculation Modes::Geometry Optimization
570 !%Description
571 !% The Fire algorithm (<tt>GOMethod = fire</tt>) uses a molecular dynamics
572 !% integrator to compute new geometries and velocities.
573 !% Currently, two integrator schemes can be selected
574 !%Option verlet 1
575 !% The Velocity Verlet algorithm.
576 !%Option euler 0
577 !% The Euler method.
578 !%End
579 call parse_variable(sys%namespace, 'GOFireIntegrator', option__gofireintegrator__verlet, g_opt%fire_integrator)
580
581 call messages_obsolete_variable(sys%namespace, 'GOWhat2Minimize', 'GOObjective')
582
583 !%Variable GOObjective
584 !%Type integer
585 !%Default minimize_energy
586 !%Section Calculation Modes::Geometry Optimization
587 !%Description
588 !% This rather esoteric option allows one to choose which
589 !% objective function to minimize during a geometry
590 !% minimization. The use of this variable may lead to
591 !% inconsistencies, so please make sure you know what you are
592 !% doing.
593 !%Option minimize_energy 1
594 !% Use the total energy as objective function.
595 !%Option minimize_forces 2
596 !% Use <math>\sqrt{\sum_i \left| f_i \right|^2}</math> as objective function.
597 !% Note that in this case one still uses the forces as the gradient of the objective function.
598 !% This is, of course, inconsistent, and may lead to very strange behavior.
599 !%End
600 call parse_variable(sys%namespace, 'GOObjective', minwhat_energy, g_opt%what2minimize)
601 if (.not. varinfo_valid_option('GOObjective', g_opt%what2minimize)) call messages_input_error(sys%namespace, 'GOObjective')
602 call messages_print_var_option("GOObjective", g_opt%what2minimize, namespace=sys%namespace)
603
604
605 !%Variable XYZGOConstrains
606 !%Type string
607 !%Section Calculation Modes::Geometry Optimization
608 !%Description
609 !% <tt>Octopus</tt> will try to read the coordinate-dependent constrains from the XYZ file
610 !% specified by the variable <tt>XYZGOConstrains</tt>.
611 !% Note: It is important for the contrains to maintain the ordering
612 !% in which the atoms were defined in the coordinates specifications.
613 !% Moreover, constrains impose fixed absolute coordinates, therefore
614 !% constrains are not compatible with GOCenter = yes
615 !%End
616
617 !%Variable XSFGOConstrains
618 !%Type string
619 !%Section Calculation Modes::Geometry Optimization
620 !%Description
621 !% Like <tt>XYZGOConstrains</tt> but in XCrySDen format, as in <tt>XSFCoordinates</tt>.
622 !%End
623
624 !%Variable PDBGOConstrains
625 !%Type string
626 !%Section Calculation Modes::Geometry Optimization
627 !%Description
628 !% Like <tt>XYZGOConstrains</tt> but in PDB format, as in <tt>PDBCoordinates</tt>.
629 !%End
630
631 !%Variable GOConstrains
632 !%Type block
633 !%Section Calculation Modes::Geometry Optimization
634 !%Description
635 !% If <tt>XYZGOConstrains</tt>, <tt>PDBConstrains</tt>, and <tt>XSFGOConstrains</tt>
636 !% are not present, <tt>Octopus</tt> will try to fetch the geometry optimization
637 !% contrains from this block. If this block is not present, <tt>Octopus</tt>
638 !% will not set any constrains. The format of this block can be
639 !% illustrated by this example:
640 !%
641 !% <tt>%GOConstrains
642 !% <br>&nbsp;&nbsp;'C' | 1 | 0 | 0
643 !% <br>&nbsp;&nbsp;'O' | &nbsp;1 | 0 | 0
644 !% <br>%</tt>
645 !%
646 !% Coordinates with a constrain value of 0 will be optimized, while
647 !% coordinates with a constrain different from zero will be kept fixed. So,
648 !% in this example the x coordinates of both atoms will remain fixed and the
649 !% distance between the two atoms along the x axis will be constant.
650 !%
651 !% Note: It is important for the constrains to maintain the ordering
652 !% in which the atoms were defined in the coordinates specifications.
653 !% Moreover, constrains impose fixed absolute coordinates, therefore
654 !% constrains are not compatible with GOCenter = yes
655 !%End
656
657 call read_coords_init(xyz)
658 call read_coords_read('GOConstrains', xyz, g_opt%ions%space, sys%namespace)
659 if (xyz%source /= read_coords_err) then
660 !Sanity check
661 if (g_opt%ions%natoms /= xyz%n) then
662 write(message(1), '(a,i4,a,i4)') 'I need exactly ', g_opt%ions%natoms, ' constrains, but I found ', xyz%n
663 call messages_fatal(1, namespace=sys%namespace)
664 end if
665 ! copy information and adjust units
666 do iatom = 1, g_opt%ions%natoms
667 where(abs(xyz%atom(iatom)%x) <= m_epsilon)
668 g_opt%ions%atom(iatom)%c = m_zero
669 elsewhere
670 g_opt%ions%atom(iatom)%c = m_one
671 end where
672 end do
673
674 call read_coords_end(xyz)
675
676 if (g_opt%fixed_atom > 0) then
677 call messages_not_implemented("GOCenter with constrains", namespace=sys%namespace)
678 end if
679 else
680 do iatom = 1, g_opt%ions%natoms
681 g_opt%ions%atom(iatom)%c = m_zero
682 end do
683 end if
684
685 call io_rm("geom/optimization.log", sys%namespace)
686
687 call io_rm("work-geom.xyz", sys%namespace)
688
689 if (.not. fromscratch) then
690 inquire(file = './last.xyz', exist = does_exist)
691 if (.not. does_exist) fromscratch = .true.
692 end if
693
694 if (.not. fromscratch) call g_opt%ions%read_xyz('./last')
695
696 ! clean out old geom/go.XXXX.xyz files. must be consistent with write_iter_info
697 iter = 1
698 do
699 write(filename, '(a,i4.4,a)') "geom/go.", iter, ".xyz"
700 inquire(file = trim(filename), exist = does_exist)
701 if (does_exist) then
702 call io_rm(trim(filename), sys%namespace)
703 iter = iter + 1
704 else
705 exit
706 end if
707 ! TODO: clean forces directory
708 end do
709
710 call restart_init(g_opt%scfv%restart_dump, sys%namespace, restart_gs, restart_type_dump, sys%mc, ierr, mesh=sys%gr)
711
713 end subroutine init_
714
715
716 ! ---------------------------------------------------------
717 subroutine end_()
718 push_sub(geom_opt_run_legacy.end_)
719
720 call states_elec_deallocate_wfns(sys%st)
721
722 call restart_end(g_opt%scfv%restart_dump)
723
724 nullify(g_opt%mesh)
725 nullify(g_opt%ions)
726 nullify(g_opt%st)
727 nullify(g_opt%hm)
728 nullify(g_opt%syst)
729
730 safe_deallocate_a(g_opt%cell_force)
731
733 end subroutine end_
734
735 end subroutine geom_opt_run_legacy
736
737
738 ! ---------------------------------------------------------
741 subroutine calc_point(size, coords, objective, getgrad, df)
742 integer, intent(in) :: size
743 real(real64), intent(in) :: coords(size)
744 real(real64), intent(inout) :: objective
745 integer, intent(in) :: getgrad
746 real(real64), intent(inout) :: df(size)
747
748 integer :: iatom, idir, jdir
749 real(real64), dimension(g_opt%periodic_dim, g_opt%periodic_dim) :: stress, scaling
750
751
752 push_sub(calc_point)
753
754 assert(size == g_opt%size)
755
756 call from_coords(g_opt, coords)
757
758 if (g_opt%fixed_atom /= 0) then
759 call g_opt%ions%translate(g_opt%ions%center())
760 end if
761
762 ! When the system is periodic in some directions, the atoms might have moved to a an adjacent cell, so we need to move them back to the original cell
763 call g_opt%ions%fold_atoms_into_cell()
764
765 ! Some atoms might have moved outside the simulation box. We stop if this happens.
766 do iatom = 1, g_opt%ions%natoms
767 if (.not. g_opt%syst%gr%box%contains_point(g_opt%syst%ions%pos(:, iatom))) then
768 if (g_opt%syst%space%periodic_dim /= g_opt%syst%space%dim) then
769 ! FIXME: This could fail for partial periodicity systems
770 ! because contains_point is too strict with atoms close to
771 ! the upper boundary to the cell.
772 write(message(1), '(a,i5,a)') "Atom ", iatom, " has moved outside the box during the geometry optimization."
773 call messages_fatal(1, namespace=g_opt%syst%namespace)
774 end if
775 end if
776 end do
777
778 call g_opt%ions%write_xyz('./work-geom', append = .true.)
779
780 call scf_mix_clear(g_opt%scfv)
781
782 ! Update lattice vectors and regenerate grid
783 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0 ) then
784 call electrons_lattice_vectors_update(g_opt%syst%namespace, g_opt%syst%gr, &
785 g_opt%syst%space, g_opt%syst%hm%psolver, g_opt%syst%hm%kpoints, &
786 g_opt%syst%mc, g_opt%syst%st%qtot, g_opt%ions%latt)
787 end if
788
789 call hamiltonian_elec_epot_generate(g_opt%hm, g_opt%syst%namespace, g_opt%syst%space, g_opt%syst%gr, &
790 g_opt%ions, g_opt%syst%ext_partners, g_opt%st)
791 call density_calc(g_opt%st, g_opt%syst%gr, g_opt%st%rho)
792 call v_ks_calc(g_opt%syst%ks, g_opt%syst%namespace, g_opt%syst%space, g_opt%hm, g_opt%st, &
793 g_opt%ions,g_opt%syst%ext_partners, calc_eigenval = .true.)
794 call energy_calc_total(g_opt%syst%namespace, g_opt%syst%space, g_opt%hm, g_opt%syst%gr, g_opt%st, g_opt%syst%ext_partners)
795
796 ! do SCF calculation
797 call scf_run(g_opt%scfv, g_opt%syst%namespace, g_opt%syst%space, g_opt%syst%mc, g_opt%syst%gr, &
798 g_opt%ions, g_opt%syst%ext_partners, &
799 g_opt%st, g_opt%syst%ks, g_opt%hm, outp = g_opt%syst%outp, verbosity = verb_compact, restart_dump=g_opt%scfv%restart_dump)
800
801 call scf_print_mem_use(g_opt%syst%namespace)
802
803 ! Convert stress into cell force
804 ! This is the Parrinello-Rahman equation of motion of the cell,
805 ! see Parrinello and Rahman, J. Appl. Pys. 52, 7182 (1981), Eq. 2.10.
806 ! Here we use a slightly different definition, in which we evolve the right stretch tensor
807 ! instead of h, because this is a symmetric matrix.
808 ! The right stretch tensor is defined as U = (1+\epsilon) = h h_0^{-1},
809 ! where \epsilon is the infinitesimal strain tensor
810 if (bitand(g_opt%type, go_cell) /= 0) then
811 stress = g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
812 do idir = 1, g_opt%periodic_dim
813 stress(idir, idir) = -stress(idir, idir) - g_opt%pressure/g_opt%periodic_dim
814 end do
815 g_opt%cell_force(:, :) = stress(:, :) * g_opt%ions%latt%rcell_volume
816
817 ! We employ here the expression 6.a of Wentzcovitch, PRB 44, 2358 (1991)
818 scaling = matmul(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim,1:g_opt%periodic_dim), g_opt%inv_initial_rlattice)
819 call lalg_inverse(g_opt%periodic_dim, scaling, 'dir')
820 g_opt%cell_force = matmul(g_opt%cell_force, transpose(scaling))
821
822 ! Thanks to the polar decomposition, any tensor can be expressed as a rotation times a symmetric matric
823 ! In order to suppress any rotation, we only keep the symmetric part
824 g_opt%cell_force = m_half * (g_opt%cell_force + transpose(g_opt%cell_force))
825 end if
826
827 ! Convert stress into cell force
828 if (bitand(g_opt%type, go_volume) /= 0) then
829 stress = g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
830 do idir = 1, g_opt%periodic_dim
831 g_opt%cell_force(idir, 1) = -(g_opt%pressure/g_opt%periodic_dim + stress(idir, idir)) * g_opt%ions%latt%rcell_volume
832 end do
833 end if
835 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0) then
836 write(message(1),'(a,3a,a)') ' Stress tensor [', trim(units_abbrev(units_out%length)), ']'
837 do idir = 1, g_opt%periodic_dim
838 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%length, g_opt%syst%st%stress_tensors%total(jdir, idir)), &
839 jdir = 1, g_opt%periodic_dim)
840 end do
841 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.true.)
842 write(message(1),'(a,3a,a)') ' Cell force tensor [', trim(units_abbrev(units_out%length)), ']'
843 do idir = 1, ubound(g_opt%cell_force, 2)
844 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%length, g_opt%cell_force(jdir, idir)), &
845 jdir = 1, g_opt%periodic_dim)
846 end do
847 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.true.)
848 end if
849
850
851 ! store results
852 if (getgrad == 1) call to_grad(g_opt, df)
853
854 if (g_opt%what2minimize == minwhat_forces) then
855 objective = m_zero
856 do iatom = 1, g_opt%ions%natoms
857 if (g_opt%ions%fixed(iatom)) cycle
858 objective = objective + sum(g_opt%ions%tot_force(:, iatom)**2)
859 end do
860 if (bitand(g_opt%type, go_cell) /= 0) then
861 do idir = 1, g_opt%periodic_dim
862 objective = objective + sum(g_opt%cell_force(:, idir)**2)
863 end do
864 end if
865 if (bitand(g_opt%type, go_volume) /= 0) then
866 objective = objective + sum(g_opt%cell_force(:,1)**2)
867 end if
868 objective = sqrt(objective)
869 else
870 objective = g_opt%hm%energy%total
871 end if
872
873 pop_sub(calc_point)
874 end subroutine calc_point
875
876
877 ! ---------------------------------------------------------
882 subroutine calc_point_ng(size, coords, objective)
883 integer :: size
884 real(real64) :: coords(size)
885 real(real64) :: objective
886
887 integer :: getgrad
888 real(real64), allocatable :: df(:)
889
890 push_sub(calc_point_ng)
891
892 assert(size == g_opt%size)
893
894 getgrad = 0
895 safe_allocate(df(1:size))
896 df = m_zero
897
898 call calc_point(size, coords, objective, getgrad, df)
899 safe_deallocate_a(df)
900
901 pop_sub(calc_point_ng)
902 end subroutine calc_point_ng
903
904
905 ! ---------------------------------------------------------
907 subroutine write_iter_info(geom_iter, size, energy, maxdx, maxdf, coords)
908 integer, intent(in) :: geom_iter
909 integer, intent(in) :: size
910 real(real64), intent(in) :: energy, maxdx, maxdf
911 real(real64), intent(in) :: coords(size)
912
913 character(len=256) :: c_geom_iter, title, c_forces_iter
914 integer :: iunit
915
916 push_sub(write_iter_info)
917
918 write(c_geom_iter, '(a,i4.4)') "go.", geom_iter
919 write(title, '(f16.10,2x,a)') units_from_atomic(units_out%energy, energy), trim(units_abbrev(units_out%energy))
920 call io_mkdir('geom', g_opt%ions%namespace)
921 call g_opt%ions%write_xyz('geom/'//trim(c_geom_iter), comment = trim(title))
922 call g_opt%ions%write_xyz('./last')
923
924 if(g_opt%periodic_dim > 0) then
925 call g_opt%ions%write_xyz('geom/'//trim(c_geom_iter), comment = 'Reduced coordinates', reduce_coordinates = .true.)
926 call write_xsf_geometry_file('geom', trim(c_geom_iter), g_opt%ions%space, g_opt%ions%latt, &
927 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace)
928 end if
929
930 if (g_opt%syst%outp%what(option__output__forces)) then
931 write(c_forces_iter, '(a,i4.4)') "forces.", geom_iter
932 if (bitand(g_opt%syst%outp%how(option__output__forces), option__outputformat__bild) /= 0) then
933 call g_opt%ions%write_bild_forces_file('forces', trim(c_forces_iter))
934 else
935 call write_xsf_geometry_file('forces', trim(c_forces_iter), g_opt%ions%space, g_opt%ions%latt, &
936 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace, total_forces=g_opt%ions%tot_force)
937 end if
938 end if
939
940 if (mpi_grp_is_root(mpi_world)) then
941 iunit = io_open(trim('geom/optimization.log'), g_opt%syst%namespace, &
942 action = 'write', position = 'append')
943
944 if (geom_iter == 1) then
945 if (bitand(g_opt%type, go_cell) /= 0) then
946 write(iunit, '(a10,5(5x,a20),a)') '# iter','energy [' // trim(units_abbrev(units_out%energy)) // ']', &
947 'max_force [' // trim(units_abbrev(units_out%force)) // ']',&
948 ' max_dr [' // trim(units_abbrev(units_out%length)) // ']', &
949 ' a, b, c ['// trim(units_abbrev(units_out%length)) // ']', &
950 ' volume ['// trim(units_abbrev(units_out%length**3)) // ']',&
951 ' alpha, beta, gamma [degrees]'
952 else
953 write(iunit, '(a10,3(5x,a20))') '# iter','energy [' // trim(units_abbrev(units_out%energy)) // ']', &
954 'max_force [' // trim(units_abbrev(units_out%force)) // ']',&
955 ' max_dr [' // trim(units_abbrev(units_out%length)) // ']'
956 end if
957 end if
958
959 if (bitand(g_opt%type, go_cell) /= 0) then
960 write(iunit, '(i10,10f25.15)') geom_iter, units_from_atomic(units_out%energy, energy), &
961 units_from_atomic(units_out%force,maxdf), &
962 units_from_atomic(units_out%length,maxdx), &
963 units_from_atomic(units_out%length,norm2(g_opt%ions%latt%rlattice(1:3, 1))),&
964 units_from_atomic(units_out%length,norm2(g_opt%ions%latt%rlattice(1:3, 2))),&
965 units_from_atomic(units_out%length,norm2(g_opt%ions%latt%rlattice(1:3, 3))),&
966 units_from_atomic(units_out%length**3, g_opt%ions%latt%rcell_volume), &
967 g_opt%ions%latt%alpha, g_opt%ions%latt%beta, g_opt%ions%latt%gamma
968 else
969 write(iunit, '(i10,3f25.15)') geom_iter, units_from_atomic(units_out%energy, energy), &
970 units_from_atomic(units_out%force,maxdf), &
971 units_from_atomic(units_out%length,maxdx)
972 end if
973
974 call io_close(iunit)
975 end if
976
977 call messages_new_line()
978 call messages_new_line()
979
980 call messages_write("++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .true.)
981
982 call messages_write("+++++++++++++++++++++ MINIMIZATION ITER #:")
983 call messages_write(geom_iter, fmt = "I5")
984 call messages_write(" ++++++++++++++++++++++", new_line = .true.)
985
986 call messages_write(" Energy = ")
987 call messages_write(energy, units = units_out%energy, fmt = "f16.10,1x", print_units = .true., new_line = .true.)
988
989 if (g_opt%periodic_dim == 0) then
990 if (maxdf > m_zero) then
991 call messages_write(" Max force = ")
992 call messages_write(maxdf, units = units_out%force, fmt = "f16.10,1x", print_units = .true., new_line = .true.)
993 end if
994
995 call messages_write(" Max dr = ")
996 call messages_write(maxdx, units = units_out%length, fmt = "f16.10,1x", print_units = .true., new_line = .true.)
997 else
998 if (maxdf > m_zero) then
999 call messages_write(" Max reduced force = ")
1000 call messages_write(maxdf, fmt = "f16.10,1x", print_units = .false., new_line = .true.)
1001 end if
1002
1003 call messages_write(" Max reduced dr = ")
1004 call messages_write(maxdx, fmt = "f16.10,1x", print_units = .false., new_line = .true.)
1005 end if
1006
1007 call messages_write("++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .true.)
1008 call messages_write("++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .true.)
1009 call messages_new_line()
1010 call messages_new_line()
1011 call messages_info()
1012
1013 pop_sub(write_iter_info)
1014 end subroutine write_iter_info
1015
1016 ! ---------------------------------------------------------
1018 subroutine to_coords(gopt, coords)
1019 type(geom_opt_t), intent(in) :: gopt
1020 real(real64), intent(out) :: coords(:)
1021
1022 integer :: iatom, idir, jdir, icoord
1023 real(real64) :: tmp_pos(gopt%dim), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1024
1025 push_sub(to_coords)
1026
1027 icoord = 1
1028 ! Ion dynamics
1029 if (bitand(g_opt%type, go_ions) /= 0) then
1030 do iatom = 1, gopt%ions%natoms
1031 if (gopt%fixed_atom == iatom) cycle
1032 if (gopt%ions%fixed(iatom)) cycle
1033 tmp_pos = gopt%ions%pos(1:gopt%dim, iatom)
1034 if (gopt%fixed_atom > 0) tmp_pos = tmp_pos - gopt%ions%pos(1:gopt%dim, gopt%fixed_atom)
1035 tmp_pos = gopt%ions%latt%cart_to_red(tmp_pos)
1036 do idir = 1, gopt%dim
1037 coords(icoord) = tmp_pos(idir)
1038 icoord = icoord + 1
1039 end do
1040 end do
1041 end if
1042
1043 ! Cell dynamics
1044 if (bitand(g_opt%type, go_cell) /= 0) then
1045 ! We compute the change in the right stretch tensor U, defined as h = U h_0
1046 strain = matmul(gopt%ions%latt%rlattice, g_opt%inv_initial_rlattice)
1047 do idir = 1, g_opt%periodic_dim
1048 do jdir = idir, g_opt%periodic_dim
1049 coords(icoord) = strain(idir, jdir)
1050 icoord = icoord + 1
1051 end do
1052 end do
1053 end if
1054
1055 ! Volume dynamics
1056 if (bitand(g_opt%type, go_volume) /= 0) then
1057 do idir = 1, g_opt%periodic_dim
1058 coords(icoord) = norm2(gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))/g_opt%initial_length(idir)
1059 icoord = icoord + 1
1060 end do
1061 end if
1062
1063
1064 pop_sub(to_coords)
1065 end subroutine to_coords
1066
1067 ! ---------------------------------------------------------
1069 subroutine to_grad(gopt, grad)
1070 type(geom_opt_t), intent(in) :: gopt
1071 real(real64), intent(out) :: grad(:)
1072
1073 integer :: iatom, idir, jdir, icoord
1074 real(real64) :: tmp_force(1:gopt%dim)
1075
1076 push_sub(to_grad)
1077
1078 icoord = 1
1079 ! Ion dynamics
1080 if (bitand(g_opt%type, go_ions) /= 0) then
1081 do iatom = 1, gopt%ions%natoms
1082 if (gopt%fixed_atom == iatom) cycle
1083 if (gopt%ions%fixed(iatom)) cycle
1084 do idir = 1, gopt%dim
1085 if (abs(gopt%ions%atom(iatom)%c(idir)) <= m_epsilon) then
1086 tmp_force(idir) = -gopt%ions%tot_force(idir, iatom)
1087 else
1088 tmp_force(idir) = m_zero
1089 end if
1090 if (gopt%fixed_atom > 0) then
1091 tmp_force(idir) = tmp_force(idir) + gopt%ions%tot_force(idir, gopt%fixed_atom)
1092 end if
1093 end do
1094 tmp_force = gopt%ions%latt%cart_to_red(tmp_force)
1095 do idir = 1, gopt%dim
1096 grad(icoord) = tmp_force(idir)
1097 icoord = icoord + 1
1098 end do
1099 end do
1100 end if
1101
1102 ! Cell dynamics
1103 if (bitand(g_opt%type, go_cell) /= 0) then
1104 do idir = 1, g_opt%periodic_dim
1105 do jdir = idir, g_opt%periodic_dim
1106 grad(icoord) = -g_opt%cell_force(idir, jdir)
1107 icoord = icoord + 1
1108 end do
1109 end do
1110 end if
1112 ! Volume dynamics
1113 if (bitand(g_opt%type, go_volume) /= 0) then
1114 do idir = 1, g_opt%periodic_dim
1115 grad(icoord) = -g_opt%cell_force(idir, 1)
1116 icoord = icoord + 1
1117 end do
1118 end if
1119
1120
1121 pop_sub(to_grad)
1122 end subroutine to_grad
1123
1124 ! ---------------------------------------------------------
1126 subroutine from_coords(gopt, coords)
1127 type(geom_opt_t), intent(inout) :: gopt
1128 real(real64), intent(in) :: coords(:)
1129
1130 integer :: iatom, idir, jdir, icoord
1131 real(real64) :: tmp_pos(gopt%dim, gopt%ions%natoms), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1132
1133 push_sub(from_coords)
1134
1135 ! Get the new reduced atomic coordinates
1136 tmp_pos = m_zero
1137 icoord = 1
1138 ! Ion dynamics
1139 if (bitand(g_opt%type, go_ions) /= 0) then
1140 do iatom = 1, gopt%ions%natoms
1141 if (gopt%fixed_atom == iatom) cycle
1142 if (gopt%ions%fixed(iatom)) cycle
1143 do idir = 1, gopt%dim
1144 tmp_pos(idir, iatom) = coords(icoord)
1145 icoord = icoord + 1
1146 end do
1147 end do
1148 else
1149 do iatom = 1, gopt%ions%natoms
1150 tmp_pos(:, iatom) = gopt%ions%latt%cart_to_red(gopt%ions%pos(:, iatom))
1151 end do
1152 end if
1153
1154 ! Updating the lattice vectors
1155 if (bitand(g_opt%type, go_cell) /= 0) then
1156 do idir = 1, g_opt%periodic_dim
1157 do jdir = idir, g_opt%periodic_dim
1158 strain(idir, jdir) = coords(icoord)
1159 icoord = icoord + 1
1160 end do
1161 end do
1162 call upper_triangular_to_hermitian(g_opt%periodic_dim, strain)
1163
1164 ! Get the new lattice vectors from the new right stretch tensor U
1165 ! The strain tensor is defined as A = U * A_0
1166 gopt%ions%latt%rlattice = matmul(strain, g_opt%initial_rlattice)
1167 end if
1168
1169 ! Updating the lattice vectors
1170 if (bitand(g_opt%type, go_volume) /= 0) then
1171 do idir = 1, g_opt%periodic_dim
1172 gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir) = coords(icoord) &
1173 * gopt%initial_rlattice(1:g_opt%periodic_dim, idir)
1174 icoord = icoord + 1
1175 end do
1176 end if
1177
1178 ! Symmetrize and update the lattice vectors
1179 if (bitand(g_opt%type, go_cell) /= 0 .or. bitand(g_opt%type, go_volume) /= 0) then
1180 call g_opt%syst%gr%symmetrizer%symmetrize_lattice_vectors(g_opt%periodic_dim, g_opt%initial_rlattice, &
1181 gopt%ions%latt%rlattice, gopt%symmetrize)
1182 call gopt%ions%update_lattice_vectors(gopt%ions%latt, gopt%symmetrize)
1183 end if
1184
1185 ! Ion dynamics
1186 if (bitand(g_opt%type, go_ions) /= 0) then
1187 ! To Cartesian coordinates
1188 do iatom = 1, gopt%ions%natoms
1189 if (gopt%fixed_atom == iatom) cycle
1190 if (gopt%ions%fixed(iatom)) cycle
1191 tmp_pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1192 do idir = 1, gopt%dim
1193 if (abs(gopt%ions%atom(iatom)%c(idir)) <= m_epsilon) then
1194 gopt%ions%pos(idir, iatom) = tmp_pos(idir, iatom)
1195 end if
1196 end do
1197 if (gopt%fixed_atom > 0) then
1198 gopt%ions%pos(:, iatom) = gopt%ions%pos(:, iatom) + gopt%ions%pos(:, gopt%fixed_atom)
1199 end if
1200 end do
1201 else
1202 do iatom = 1, gopt%ions%natoms
1203 gopt%ions%pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1204 end do
1205 end if
1206
1207 if (gopt%symmetrize) then
1208 call gopt%ions%symmetrize_atomic_coord()
1209 end if
1210
1211 if (debug%info) then
1212 call gopt%ions%print_spacegroup()
1213 end if
1214
1215 pop_sub(from_coords)
1216 end subroutine from_coords
1217
1218 ! ---------------------------------------------------------
1220 subroutine write_iter_info_ng(geom_iter, size, energy, maxdx, coords)
1221 integer, intent(in) :: geom_iter
1222 integer, intent(in) :: size
1223 real(real64), intent(in) :: energy, maxdx
1224 real(real64), intent(in) :: coords(size)
1225
1226 push_sub(write_iter_info_ng)
1227 call write_iter_info(geom_iter, size, energy, maxdx, -m_one, coords)
1228
1229 pop_sub(write_iter_info_ng)
1230 end subroutine write_iter_info_ng
1231
1232end module geom_opt_oct_m
1233
1234!! Local Variables:
1235!! mode: f90
1236!! coding: utf-8
1237!! End:
subroutine init_(fromscratch)
Definition: geom_opt.F90:354
subroutine end_()
Definition: geom_opt.F90:811
Define which routines can be seen from the outside.
Definition: loct.F90:147
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
Definition: minimizer.F90:403
subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
Implementation of the Fast Inertial Relaxation Engine (FIRE)
Definition: minimizer.F90:482
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
type(debug_t), save, public debug
Definition: debug.F90:156
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
subroutine, public geom_opt_run(system, from_scratch)
Definition: geom_opt.F90:207
subroutine calc_point_ng(size, coords, objective)
Same as calc_point, but without the gradients. No intents here is unfortunately required because the ...
Definition: geom_opt.F90:976
subroutine to_grad(gopt, grad)
Transfer data from the forces to the work array for the gradients (grad)
Definition: geom_opt.F90:1163
integer, parameter go_cell
Definition: geom_opt.F90:198
integer, parameter minwhat_forces
Definition: geom_opt.F90:194
subroutine write_iter_info_ng(geom_iter, size, energy, maxdx, coords)
Same as write_iter_info, but without the gradients.
Definition: geom_opt.F90:1314
subroutine write_iter_info(geom_iter, size, energy, maxdx, maxdf, coords)
Output the information after each iteration of the geometry optimization.
Definition: geom_opt.F90:1001
subroutine calc_point(size, coords, objective, getgrad, df)
Note: you might think it would be better to change the arguments with '(size)' below to '(:)'....
Definition: geom_opt.F90:835
subroutine to_coords(gopt, coords)
Transfer the data from the data structures to the work array (coords)
Definition: geom_opt.F90:1112
integer, parameter go_volume
Definition: geom_opt.F90:198
subroutine from_coords(gopt, coords)
Transfer the data from the work array (coords) to the actual data structures.
Definition: geom_opt.F90:1220
subroutine geom_opt_run_legacy(sys, fromscratch)
Definition: geom_opt.F90:225
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_rm(fname, namespace)
Definition: io.F90:342
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
Definition: lcao.F90:808
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1113
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
subroutine, public messages_new_line()
Definition: messages.F90:1134
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1085
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
integer minmethod_fire
Definition: minimizer.F90:134
integer minmethod_nmsimplex
Definition: minimizer.F90:134
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module implements the basic mulsisystem class, a container system for other systems.
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer, parameter, public read_coords_err
for read_coords_info::file_type
subroutine, public read_coords_init(gf)
subroutine, public read_coords_end(gf)
subroutine, public read_coords_read(what, gf, space, namespace)
integer, parameter, public restart_gs
Definition: restart.F90:200
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:516
integer, parameter, public restart_type_dump
Definition: restart.F90:245
integer, parameter, public restart_type_load
Definition: restart.F90:245
subroutine, public restart_end(restart)
Definition: restart.F90:722
subroutine, public scf_print_mem_use(namespace)
Definition: scf.F90:1575
subroutine, public scf_mix_clear(scf)
Definition: scf.F90:577
integer, parameter, public verb_compact
Definition: scf.F90:203
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
Definition: scf.F90:254
subroutine, public scf_end(scf)
Definition: scf.F90:547
subroutine, public scf_run(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, verbosity, iters_done, restart_dump)
Legacy version of the SCF code.
Definition: scf.F90:826
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_t), public unit_femtosecond
Time in femtoseconds.
type(unit_t), public unit_amu
Mass in atomic mass units (AKA Dalton).
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:736
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
Definition: v_ks.F90:683
An abstract type for all electron species.
Class describing the electron system.
Definition: electrons.F90:218
Container class for lists of system_oct_m::system_t.
int true(void)