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