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