Octopus
ion_dynamics.F90
Go to the documentation of this file.
1!! Copyright (C) 2008 X. Andrade
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
22 use accel_oct_m
24 use iso_c_binding
25 use debug_oct_m
26 use global_oct_m
27 use grid_oct_m
28 use ions_oct_m
29 use, intrinsic :: iso_fortran_env
34 use math_oct_m
36 use mpi_oct_m
39 use parser_oct_m
43 use space_oct_m
48 use unit_oct_m
51
52 implicit none
53
54 private
55
56 public :: &
75
76
77 integer, parameter :: &
78 THERMO_NONE = 0, &
79 thermo_scal = 1, &
80 thermo_nh = 2
81
82 type nose_hoover_t
83 private
84 real(real64) :: mass
85 real(real64) :: pos
86 real(real64) :: vel
87 end type nose_hoover_t
88
90 private
91 logical :: move
92 type(tdf_t) :: fx
93 type(tdf_t) :: fy
94 type(tdf_t) :: fz
96
98 private
99 logical :: move_ions
100 logical :: relax_cell
101 logical :: constant_velocity
102 integer :: thermostat
103 real(real64) :: dt
104 real(real64) :: current_temperature
105
106 real(real64), allocatable :: oldforce(:, :)
107
109 real(real64), allocatable :: old_pos(:, :)
110
112 real(real64), allocatable :: cell_force(:)
113 real(real64), allocatable :: old_cell_force(:)
114 real(real64), allocatable :: cell_vel(:)
115 real(real64), allocatable :: initial_rlattice(:,:)
116 real(real64), allocatable :: strain(:)
117
118 real(real64) :: pressure
119
121 logical :: symmetrize = .false.
122 type(symmetrizer_t), pointer :: symm
123
125 type(nose_hoover_t) :: nh(1:2)
126 type(tdf_t) :: temperature_function
127
129 logical :: drive_ions
130 type(ion_td_displacement_t), allocatable :: td_displacements(:)
131 type(ions_t), pointer :: ions_t0
132
133 real(real64), public :: ionic_scale
134 contains
135 procedure :: ions_move => ion_dynamics_ions_move
136 procedure :: cell_relax => ion_dynamics_cell_relax
137 procedure :: update_stress => ion_dynamics_update_stress
138 procedure :: is_active => ion_dynamics_is_active
139 end type ion_dynamics_t
140
141 type ion_state_t
142 private
143 real(real64), allocatable :: pos(:, :)
144 real(real64), allocatable :: vel(:, :)
145 real(real64), allocatable :: old_pos(:, :)
146 type(nose_hoover_t) :: nh(1:2)
147 end type ion_state_t
148
149 type cell_state_t
150 private
151 real(real64), allocatable :: pos(:, :)
152 real(real64), allocatable :: vel(:, :)
153 real(real64), allocatable :: old_pos(:, :)
154 end type cell_state_t
155
156contains
157
158 ! ---------------------------------------------------------
159 subroutine ion_dynamics_init(this, namespace, ions, symmetrize, symm)
160 use iso_c_binding, only: c_ptr
161 type(ion_dynamics_t), intent(out) :: this
162 type(namespace_t), intent(in) :: namespace
163 type(ions_t), intent(inout) :: ions
164 logical, intent(in) :: symmetrize
165 type(symmetrizer_t), optional, target, intent(in) :: symm
166
167 integer :: i, j, iatom, ierr, periodic_dim, ncomp
168 real(real64) :: xx(ions%space%dim), temperature, sigma, kin1, kin2
169 type(c_ptr) :: random_gen_pointer
170 type(read_coords_info) :: xyz
171 character(len=100) :: temp_function_name
172 logical :: have_velocities
173
174 type(block_t) :: blk
175 integer :: ndisp
176 character(len=200) :: expression
180 have_velocities = .false.
181 this%drive_ions = .false.
183 this%symmetrize = symmetrize
184 if (this%symmetrize) then
185 assert(present(symm))
186 this%symm => symm
187 else
188 nullify(this%symm)
189 end if
191 !%Variable TDIonicTimeScale
192 !%Type float
193 !%Default 1.0
194 !%Section Time-Dependent::Propagation
195 !%Description
196 !% This variable defines the factor between the timescale of ionic
197 !% and electronic movement. It allows reasonably fast
198 !% Born-Oppenheimer molecular-dynamics simulations based on
199 !% Ehrenfest dynamics. The value of this variable is equivalent to
200 !% the role of <math>\mu</math> in Car-Parrinello. Increasing it
201 !% linearly accelerates the time step of the ion
202 !% dynamics, but also increases the deviation of the system from the
203 !% Born-Oppenheimer surface. The default is 1, which means that both
204 !% timescales are the same. Note that a value different than 1
205 !% implies that the electrons will not follow physical behaviour.
206 !%
207 !% According to our tests, values around 10 are reasonable, but it
208 !% will depend on your system, mainly on the width of the gap.
209 !%
210 !% Important: The electronic time step will be the value of
211 !% <tt>TDTimeStep</tt> divided by this variable, so if you have determined an
212 !% optimal electronic time step (that we can call <i>dte</i>), it is
213 !% recommended that you define your time step as:
214 !%
215 !% <tt>TDTimeStep</tt> = <i>dte</i> * <tt>TDIonicTimeScale</tt>
216 !%
217 !% so you will always use the optimal electronic time step
218 !% (<a href=http://arxiv.org/abs/0710.3321>more details</a>).
219 !%End
220 call parse_variable(namespace, 'TDIonicTimeScale', m_one, this%ionic_scale)
221
222 if (this%ionic_scale <= m_zero) then
223 write(message(1),'(a)') 'Input: TDIonicTimeScale must be positive.'
224 call messages_fatal(1, namespace=namespace)
225 end if
227 call messages_print_var_value('TDIonicTimeScale', this%ionic_scale, namespace=namespace)
232 !%Variable IonsConstantVelocity
233 !%Type logical
234 !%Default no
235 !%Section Time-Dependent::Propagation
236 !%Description
237 !% (Experimental) If this variable is set to yes, the ions will
238 !% move with a constant velocity given by the initial
239 !% conditions. They will not be affected by any forces.
240 !%End
241 call parse_variable(namespace, 'IonsConstantVelocity', .false., this%constant_velocity)
242 call messages_print_var_value('IonsConstantVelocity', this%constant_velocity, namespace=namespace)
243
244 if (this%constant_velocity) then
245 call messages_experimental('IonsConstantVelocity', namespace=namespace)
246 have_velocities = .true.
247 this%drive_ions = .true.
248 end if
249
250 !%Variable IonsTimeDependentDisplacements
251 !%Type block
252 !%Section Time-Dependent::Propagation
253 !%Description
254 !% (Experimental) This variable allows you to specify a
255 !% time-dependent function describing the displacement of the ions
256 !% from their equilibrium position: <math>r(t) = r_0 + \Delta
257 !% r(t)</math>. Specify the displacements dx(t), dy(t), dz(t) as
258 !% follows, for some or all of the atoms:
259 !%
260 !% <tt>%IonsTimeDependentDisplacements
261 !% <br>&nbsp;&nbsp; atom_index | "dx(t)" | "dy(t)" | "dz(t)"
262 !% <br>%</tt>
263 !%
264 !% The displacement functions are time-dependent functions and should match one
265 !% of the function names given in the first column of the <tt>TDFunctions</tt> block.
266 !% If this block is set, the ions will not be affected by any forces.
267 !%End
268
269
270 ndisp = 0
271 if (parse_block(namespace, 'IonsTimeDependentDisplacements', blk) == 0) then
272 call messages_experimental("IonsTimeDependentDisplacements", namespace=namespace)
273 ndisp= parse_block_n(blk)
274 safe_allocate(this%td_displacements(1:ions%natoms))
275 this%td_displacements(1:ions%natoms)%move = .false.
276 if (ndisp > 0) this%drive_ions =.true.
277
278 do i = 1, ndisp
279 call parse_block_integer(blk, i-1, 0, iatom)
280 this%td_displacements(iatom)%move = .true.
281
282 call parse_block_string(blk, i-1, 1, expression)
283 call tdf_read(this%td_displacements(iatom)%fx, namespace, trim(expression), ierr)
284 if (ierr /= 0) then
285 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
286 call messages_warning(1, namespace=namespace)
287 end if
288
289
290 call parse_block_string(blk, i-1, 2, expression)
291 call tdf_read(this%td_displacements(iatom)%fy, namespace, trim(expression), ierr)
292 if (ierr /= 0) then
293 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
294 call messages_warning(1, namespace=namespace)
295 end if
296
297 call parse_block_string(blk, i-1, 3, expression)
298 call tdf_read(this%td_displacements(iatom)%fz, namespace, trim(expression), ierr)
299 if (ierr /= 0) then
300 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
301 call messages_warning(1, namespace=namespace)
302 end if
303
304 end do
305
306 safe_allocate(this%ions_t0)
307 this%ions_t0 = ions
308
309 end if
310
311
312
313
314 !%Variable Thermostat
315 !%Type integer
316 !%Default none
317 !%Section Time-Dependent::Propagation
318 !%Description
319 !% This variable selects the type of thermostat applied to
320 !% control the ionic temperature.
321 !%Option none 0
322 !% No thermostat is applied. This is the default.
323 !%Option velocity_scaling 1
324 !% Velocities are scaled to control the temperature.
325 !%Option nose_hoover 2
326 !% Nose-Hoover thermostat.
327 !%End
328
329 call parse_variable(namespace, 'Thermostat', thermo_none, this%thermostat)
330 if (.not. varinfo_valid_option('Thermostat', this%thermostat)) call messages_input_error(namespace, 'Thermostat')
331 call messages_print_var_option('Thermostat', this%thermostat, namespace=namespace)
332
333 if (this%thermostat /= thermo_none) then
334
335 have_velocities = .true.
336
337 if (this%drive_ions) then
338 call messages_write('You cannot use a Thermostat and IonsConstantVelocity or IonsTimeDependentDisplacements')
339 call messages_write('at the same time.')
340 call messages_fatal(namespace=namespace)
341 end if
342
343 call messages_experimental('Thermostat', namespace=namespace)
344
345 !%Variable TemperatureFunction
346 !%Type integer
347 !%Default "temperature"
348 !%Section Time-Dependent::Propagation
349 !%Description
350 !% If a thermostat is used, this variable indicates the name of the
351 !% function in the <tt>TDFunctions</tt> block that will be used to control the
352 !% temperature. The values of the temperature are given in
353 !% degrees Kelvin.
354 !%End
355 call parse_variable(namespace, 'TemperatureFunction', 'temperature', temp_function_name)
356
357 call tdf_read(this%temperature_function, namespace, temp_function_name, ierr)
358
359 if (ierr /= 0) then
360 message(1) = "You have enabled a thermostat but Octopus could not find"
361 message(2) = "the '"//trim(temp_function_name)//"' function in the TDFunctions block."
362 call messages_fatal(2, namespace=namespace)
363 end if
364
365 if (this%thermostat == thermo_nh) then
366 !%Variable ThermostatMass
367 !%Type float
368 !%Default 1.0
369 !%Section Time-Dependent::Propagation
370 !%Description
371 !% This variable sets the fictitious mass for the Nose-Hoover
372 !% thermostat.
373 !%End
374 call messages_obsolete_variable(namespace, 'NHMass', 'ThermostatMass')
375
376 call parse_variable(namespace, 'ThermostatMass', m_one, this%nh(1)%mass)
377 this%nh(2)%mass = this%nh(1)%mass
378
379 this%nh(1:2)%pos = m_zero
380 this%nh(1:2)%vel = m_zero
381
382 safe_allocate(this%old_pos(1:ions%space%dim, 1:ions%natoms))
383
384 this%old_pos = ions%pos
385 end if
386
387 end if
388
389 !now initialize velocities
390
391 !%Variable RandomVelocityTemp
392 !%Type float
393 !%Default 0.0
394 !%Section System::Velocities
395 !%Description
396 !% If this variable is present, <tt>Octopus</tt> will assign random
397 !% velocities to the atoms following a Boltzmann distribution with
398 !% temperature given by <tt>RandomVelocityTemp</tt> (in degrees Kelvin).
399 !% The seed for the random number generator can be modified by setting
400 !% <tt>GSL_RNG_SEED</tt> environment variable.
401 !%End
402
403 ! we now load the velocities, either from the temperature, from the input, or from a file
404 if (parse_is_defined(namespace, 'RandomVelocityTemp')) then
405
406 have_velocities = .true.
407
408 if (mpi_grp_is_root(mpi_world)) then
409 call loct_ran_init(random_gen_pointer)
410 call parse_variable(namespace, 'RandomVelocityTemp', m_zero, temperature, unit = unit_kelvin)
411 end if
412
413 do i = 1, ions%natoms
414 !generate the velocities in the root node
415 if (mpi_grp_is_root(mpi_world)) then
416 sigma = sqrt(temperature / ions%mass(i))
417 do j = 1, 3
418 ions%vel(j, i) = loct_ran_gaussian(random_gen_pointer, sigma)
419 end do
420 end if
421 !and send them to the others
422 call mpi_world%bcast(ions%vel(:, i), ions%space%dim, mpi_double_precision, 0)
423 end do
424
425 if (mpi_grp_is_root(mpi_world)) then
426 call loct_ran_end(random_gen_pointer)
427 end if
428
429 call ions%update_kinetic_energy()
430 kin1 = ions%kinetic_energy
431
432 xx = ions%center_of_mass_vel()
433 do i = 1, ions%natoms
434 ions%vel(:, i) = ions%vel(:, i) - xx
435 end do
436
437 call ions%update_kinetic_energy()
438 kin2 = ions%kinetic_energy
439
440 do i = 1, ions%natoms
441 ions%vel(:, i) = sqrt(kin1/kin2)*ions%vel(:, i)
442 end do
443
444 call ions%update_kinetic_energy()
445
446 write(message(1),'(a,f10.4,1x,a)') 'Info: Initial velocities randomly distributed with T =', &
448 write(message(2),'(2x,a,f8.4,1x,a)') '<K> =', &
449 units_from_atomic(units_out%energy, ions%kinetic_energy/ions%natoms), &
450 units_abbrev(units_out%energy)
451 write(message(3),'(2x,a,f8.4,1x,a)') '3/2 k_B T =', &
452 units_from_atomic(units_out%energy, (m_three/m_two)*temperature), &
453 units_abbrev(units_out%energy)
454 call messages_info(3, namespace=namespace)
455
456 else
457 !%Variable XYZVelocities
458 !%Type string
459 !%Section System::Velocities
460 !%Description
461 !% <tt>Octopus</tt> will try to read the starting velocities of the atoms from the XYZ file
462 !% specified by the variable <tt>XYZVelocities</tt>.
463 !% Note that you do not need to specify initial velocities if you are not going
464 !% to perform ion dynamics; if you are going to allow the ions to move but the velocities
465 !% are not specified, they are considered to be null.
466 !% Note: It is important for the velocities to maintain the ordering
467 !% in which the atoms were defined in the coordinates specifications.
468 !%End
469
470 !%Variable XSFVelocities
471 !%Type string
472 !%Section System::Velocities
473 !%Description
474 !% Like <tt>XYZVelocities</tt> but in XCrySDen format, as in <tt>XSFCoordinates</tt>.
475 !%End
476
477 !%Variable PDBVelocities
478 !%Type string
479 !%Section System::Velocities
480 !%Description
481 !% Like <tt>XYZVelocities</tt> but in PDB format, as in <tt>PDBCoordinates</tt>.
482 !%End
483
484 !%Variable Velocities
485 !%Type block
486 !%Section System::Velocities
487 !%Description
488 !% If <tt>XYZVelocities</tt>, <tt>PDBVelocities</tt>, and <tt>XSFVelocities</tt>
489 !% are not present, <tt>Octopus</tt> will try to fetch the initial
490 !% atomic velocities from this block. If this block is not present, <tt>Octopus</tt>
491 !% will set the initial velocities to zero. The format of this block can be
492 !% illustrated by this example:
493 !%
494 !% <tt>%Velocities
495 !% <br>&nbsp;&nbsp;'C' | -1.7 | 0.0 | 0.0
496 !% <br>&nbsp;&nbsp;'O' | &nbsp;1.7 | 0.0 | 0.0
497 !% <br>%</tt>
498 !%
499 !% It describes one carbon and one oxygen moving at the relative
500 !% velocity of 3.4 velocity units.
501 !%
502 !% Note: It is important for the velocities to maintain the ordering
503 !% in which the atoms were defined in the coordinates specifications.
504 !%End
505
506 call read_coords_init(xyz)
507 call read_coords_read('Velocities', xyz, ions%space, namespace)
508 if (xyz%source /= read_coords_err) then
509
510 have_velocities = .true.
511
512 if (ions%natoms /= xyz%n) then
513 write(message(1), '(a,i4,a,i4)') 'I need exactly ', ions%natoms, ' velocities, but I found ', xyz%n
514 call messages_fatal(1, namespace=namespace)
515 end if
516
517 ! copy information and adjust units
518 do i = 1, ions%natoms
519 ions%vel(:, i) = units_to_atomic(units_inp%velocity/units_inp%length, xyz%atom(i)%x(1:ions%space%dim))
520 end do
521
522 call read_coords_end(xyz)
523
524 else
525 ions%vel = m_zero
526 end if
527 end if
528
529 call ions%update_kinetic_energy()
530
531 !%Variable MoveIons
532 !%Type logical
533 !%Section Time-Dependent::Propagation
534 !%Description
535 !% This variable controls whether atoms are moved during a time
536 !% propagation run. The default is yes when the ion velocity is
537 !% set explicitly or implicitly, otherwise is no.
538 !%End
539 call parse_variable(namespace, 'MoveIons', have_velocities, this%move_ions)
540 call messages_print_var_value('MoveIons', this%move_ions, namespace=namespace)
541
542 if (this%move_ions .and. ions%space%periodic_dim == 1) then
543 call messages_input_error(namespace, 'MoveIons', &
544 'Moving ions for a 1D periodic system is not allowed, as forces are incorrect.')
545 end if
546
547 if (this%ions_move()) then
548 safe_allocate(this%oldforce(1:ions%space%dim, 1:ions%natoms))
549 end if
550
551 if (ions%space%is_periodic()) then
552 !%Variable CellDynamics
553 !%Type logical
554 !%Section Time-Dependent::Propagation
555 !%Description
556 !% This variable controls whether the cell relaxation is done during a time
557 !% propagation run. The default is no.
558 !% This is done based on the Parrinello-Rahman equation of motion of the cell,
559 !% see Parrinello and Rahman, J. Appl. Pys. 52, 7182 (1981).
560 !% Only for periodic systems.
561 !%End
562 call parse_variable(namespace, 'CellDynamics', .false., this%relax_cell)
563 call messages_print_var_value('CellDynamics', this%relax_cell, namespace=namespace)
564
565 if (this%cell_relax()) then
566 if (accel_is_enabled()) then
567 message(1) = "Cell dynamics not supported on GPUs."
568 call messages_fatal(1, namespace=namespace)
569 end if
570
571 periodic_dim = ions%space%periodic_dim
572 ncomp = periodic_dim * periodic_dim
573 safe_allocate(this%cell_force(1:ncomp))
574 this%cell_force = m_zero
575 safe_allocate(this%old_cell_force(1:ncomp))
576 this%old_cell_force = m_zero
577 safe_allocate(this%cell_vel(1:ncomp))
578 this%cell_vel = m_zero
579 ! We start from identity for the strain
580 safe_allocate(this%strain(1:ncomp))
581 this%strain = m_zero
582 ncomp = 1
583 do i = 1, periodic_dim
584 do j = i, periodic_dim
585 if(i == j) this%strain(ncomp) = m_one
586 ncomp = ncomp + 1
587 end do
588 end do
589
590 ! As we work with "reduced" quantities, the initial lattice vectors are used as a reference
591 ! and the strain is propagated from this initial reference
592 safe_allocate(this%initial_rlattice(1:periodic_dim, 1:periodic_dim))
593 this%initial_rlattice(1:periodic_dim, 1:periodic_dim) = ions%latt%rlattice(1:periodic_dim, 1:periodic_dim)
594 end if
595
596 !%Variable HydrostaticPressure
597 !%Type float
598 !%Default 0.0
599 !%Section Time-Dependent::Propagation
600 !%Description
601 !% Geometry optimization and molecular dynamics can be performed in presence of an external
602 !% hydrostatic pressure. This variable allows to set this value.
603 !% Only for periodic systems.
604 !%End
605 call parse_variable(namespace, 'HydrostaticPressure', m_zero, this%pressure)
606 else
607 this%relax_cell = .false.
608 end if
609
610
611 pop_sub(ion_dynamics_init)
612 end subroutine ion_dynamics_init
613
614
615 ! ---------------------------------------------------------
616 subroutine ion_dynamics_end(this)
617 type(ion_dynamics_t), intent(inout) :: this
618
619 push_sub(ion_dynamics_end)
620 safe_deallocate_a(this%oldforce)
621
622 if (this%thermostat /= thermo_none) then
623 call tdf_end(this%temperature_function)
624 end if
625
626 if (this%drive_ions .and. allocated(this%td_displacements)) then
627 if (any(this%td_displacements(1:this%ions_t0%natoms)%move)) then
628 ! ions end cannot be called here, otherwise the species are destroyed twice
629 ! call ions_end(this%ions_t0)
630 !AFE_DEALLOCATE_P(this%ions_t0)
631 end if
632 safe_deallocate_a(this%td_displacements)
633 end if
634
635 safe_deallocate_a(this%cell_force)
636 safe_deallocate_a(this%old_cell_force)
637 safe_deallocate_a(this%cell_vel)
638 safe_deallocate_a(this%initial_rlattice)
639
640 pop_sub(ion_dynamics_end)
641 end subroutine ion_dynamics_end
642
643
644 ! ---------------------------------------------------------
646 subroutine ion_dynamics_propagate(this, ions, time, dt, namespace)
647 type(ion_dynamics_t), intent(inout) :: this
648 type(ions_t), intent(inout) :: ions
649 real(real64), intent(in) :: time
650 real(real64), intent(in) :: dt
651 type(namespace_t), intent(in) :: namespace
652
653 integer :: iatom
654
655 push_sub(ion_dynamics_propagate)
656
657 this%dt = dt
658
659 if (this%drive_ions) then
660
661 call ion_dynamics_propagate_driven_ions(this, ions, time, dt)
662
663 else
664
665 ! Update the thermostat temperature
666 call ion_dynamics_update_temperature(this, time, namespace)
667
668 ! Conversion to reduced coordinates
669 do iatom = 1, ions%natoms
670 ions%pos(:, iatom) = ions%latt%cart_to_red(ions%pos(:, iatom))
671 ions%vel(:, iatom) = ions%latt%cart_to_red(ions%vel(:, iatom))
672 ions%tot_force(:, iatom) = ions%latt%cart_to_red(ions%tot_force(:, iatom))
673 end do
674
675 ! Get the new reduced coordinates
676 if (this%ions_move()) then
677 call ion_dynamics_propagate_ions(this, ions, time, dt)
678 end if
679
680 ! Updating the lattice vectors
681 if (this%cell_relax()) then
682 call ion_dynamics_propagate_cell(this, ions, time, dt, namespace)
683 end if
684
685 ! Get the new Cartesian coordinates
686 do iatom = 1, ions%natoms
687 ions%pos(:, iatom) = ions%latt%red_to_cart(ions%pos(:, iatom))
688 ions%vel(:, iatom) = ions%latt%red_to_cart(ions%vel(:, iatom))
689 if (allocated(this%oldforce)) then
690 this%oldforce(:, iatom) = ions%latt%red_to_cart(this%oldforce(:, iatom))
691 end if
692 ions%tot_force(:, iatom) = ions%latt%red_to_cart(ions%tot_force(:, iatom))
693 end do
694
695 end if
696
697 call ions%fold_atoms_into_cell()
698
700 end subroutine ion_dynamics_propagate
701
702
703 ! ---------------------------------------------------------
705 subroutine ion_dynamics_update_temperature(this, time, namespace)
706 type(ion_dynamics_t), intent(inout) :: this
707 real(real64), intent(in) :: time
708 type(namespace_t), intent(in) :: namespace
711
712 ! get the temperature from the tdfunction for the current time
713 if (this%thermostat /= thermo_none) then
714 this%current_temperature = units_to_atomic(unit_kelvin, tdf(this%temperature_function, time))
715
716 if (this%current_temperature < m_zero) then
717 write(message(1), '(a, f10.3, 3a, f10.3, 3a)') &
718 "Negative temperature (", &
719 units_from_atomic(unit_kelvin, this%current_temperature), " ", units_abbrev(unit_kelvin), &
720 ") at time ", &
721 units_from_atomic(units_out%time, time), " ", trim(units_abbrev(units_out%time)), "."
722 call messages_fatal(1, namespace=namespace)
723 end if
724 else
725 this%current_temperature = m_zero
726 end if
727
730
734 subroutine ion_dynamics_propagate_driven_ions(this, ions, time, dt)
735 type(ion_dynamics_t), intent(inout) :: this
736 type(ions_t), intent(inout) :: ions
737 real(real64), intent(in) :: time
738 real(real64), intent(in) :: dt
740 integer :: iatom
741 real(real64) :: dr(3)
742
744
745 assert(this%drive_ions)
746
747 do iatom = 1, ions%natoms
748 if (ions%fixed(iatom)) cycle
749
750 if (this%constant_velocity) then
751 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom)
752
753 else if (allocated(this%td_displacements)) then
754
755 if (this%td_displacements(iatom)%move) then
756 dr(1:3)=(/ real(tdf(this%td_displacements(iatom)%fx, time), real64), &
757 real(tdf(this%td_displacements(iatom)%fy, time), real64), &
758 real(tdf(this%td_displacements(iatom)%fz, time), real64) /)
759
760 ions%pos(:, iatom) = this%ions_t0%pos(:, iatom) + dr(1:ions%space%dim)
761 end if
762
763 end if
764 end do
765
768
769 ! ---------------------------------------------------------
774 subroutine ion_dynamics_propagate_ions(this, ions, time, dt)
775 type(ion_dynamics_t), intent(inout) :: this
776 type(ions_t), intent(inout) :: ions
777 real(real64), intent(in) :: time
778 real(real64), intent(in) :: dt
779
780 integer :: iatom
781
783
784 assert(.not. this%drive_ions)
785
786 if (this%thermostat /= thermo_nh) then
787 ! integrate using verlet
788 do iatom = 1, ions%natoms
789 if (ions%fixed(iatom)) cycle
790
791 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom) + &
792 m_half*dt**2 / ions%mass(iatom) * ions%tot_force(:, iatom)
793
794 this%oldforce(:, iatom) = ions%tot_force(:, iatom)
795 end do
796
797 else
798 ! for the Nose-Hoover thermostat we use a special integrator
799
800 ! The implementation of the Nose-Hoover thermostat is based on
801 ! Understanding Molecular Simulations by Frenkel and Smit,
802 ! Appendix E, page 540-542.
803
804 call nh_chain(this, ions)
805
806 do iatom = 1, ions%natoms
807 if (ions%fixed(iatom)) cycle
808
809 ions%pos(:, iatom) = ions%pos(:, iatom) + m_half*dt*ions%vel(:, iatom)
810 end do
811
812 end if
813
815 end subroutine ion_dynamics_propagate_ions
816
817 ! ---------------------------------------------------------
819 subroutine ion_dynamics_propagate_cell(this, ions, time, dt, namespace)
820 type(ion_dynamics_t), intent(inout) :: this
821 type(ions_t), intent(inout) :: ions
822 real(real64), intent(in) :: time
823 real(real64), intent(in) :: dt
824 type(namespace_t), intent(in) :: namespace
825
826 integer :: idir, jdir, comp
827 real(real64) :: rlattice_change(ions%space%periodic_dim*ions%space%periodic_dim)
828
830
831 rlattice_change = dt * this%cell_vel + m_half*dt**2 * this%cell_force
832
833 write(message(1),'(a,3a,a)') ' Cell force [', trim(units_abbrev(units_out%energy/units_out%length**ions%space%dim)), ']'
834 do idir = 1, ions%space%periodic_dim
835 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%energy/units_out%length**ions%space%dim, &
836 this%cell_force(jdir + (idir-1)*ions%space%periodic_dim)), &
837 jdir = 1, ions%space%periodic_dim)
838 end do
839 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
840
841 write(message(1),'(a,3a,a)') ' Cell vel [', &
842 trim(units_abbrev(units_out%energy/units_out%length**ions%space%dim*units_out%time)), ']'
843 do idir = 1, ions%space%periodic_dim
844 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%energy/units_out%length**ions%space%dim*units_out%time, &
845 this%cell_vel(jdir+ (idir-1)*ions%space%periodic_dim)), &
846 jdir = 1, ions%space%periodic_dim)
847 end do
848 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
849
850
851 comp = 1
852 do idir = 1, ions%space%periodic_dim
853 do jdir = 1, ions%space%periodic_dim
854 ions%latt%rlattice(idir, jdir) = ions%latt%rlattice(idir, jdir) + rlattice_change(comp)
855 comp = comp + 1
856 end do
857 end do
858
859 this%old_cell_force = this%cell_force
860
861 if (associated(this%symm)) then
862 call this%symm%symmetrize_lattice_vectors(ions%space%periodic_dim, &
863 this%initial_rlattice, ions%latt%rlattice, this%symmetrize)
864 end if
865 call ions%update_lattice_vectors(ions%latt, this%symmetrize)
866
868 end subroutine ion_dynamics_propagate_cell
869
870
871 ! ---------------------------------------------------------
872 subroutine nh_chain(this, ions)
873 type(ion_dynamics_t), intent(inout) :: this
874 type(ions_t), intent(inout) :: ions
875
876 real(real64) :: g1, g2, ss, uk, dt, temp
877
878 push_sub(nh_chain)
879
880 dt = this%dt
881
882 call ions%update_kinetic_energy()
883 uk = ions%kinetic_energy
884
885 temp = this%current_temperature
886
887 g2 = (this%nh(1)%mass*this%nh(1)%vel**2 - temp)/this%nh(2)%mass
888 this%nh(2)%vel = this%nh(2)%vel + g2*dt/m_four
889 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
890
891 g1 = (m_two*uk - m_three*ions%natoms*temp)/this%nh(1)%mass
892 this%nh(1)%vel = this%nh(1)%vel + g1*dt/m_four
893 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
894 this%nh(1)%pos = this%nh(1)%pos + this%nh(1)%vel*dt/m_two
895 this%nh(2)%pos = this%nh(2)%pos + this%nh(2)%vel*dt/m_two
896
897 ss = exp(-this%nh(1)%vel*dt/m_two)
898
899 ions%vel = ss*ions%vel
900
901 uk = uk*ss**2
902
903 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
904 g1 = (m_two*uk - m_three*ions%natoms*temp)/this%nh(1)%mass
905 this%nh(1)%vel = this%nh(1)%vel + g1*dt/m_four
906 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
907
908 g2 = (this%nh(1)%mass*this%nh(1)%vel**2 - temp)/this%nh(2)%mass
909 this%nh(2)%vel = this%nh(2)%vel + g2*dt/m_four
910
911 pop_sub(nh_chain)
912 end subroutine nh_chain
913
914
915 ! ---------------------------------------------------------
916 subroutine ion_dynamics_propagate_vel(this, ions, atoms_moved)
917 type(ion_dynamics_t), intent(inout) :: this
918 type(ions_t), intent(inout) :: ions
919 logical, optional, intent(out) :: atoms_moved
920
921 integer :: iatom
922 real(real64) :: scal
923
924 if (.not. ion_dynamics_ions_move(this)) return
925 if (this%drive_ions) return
926
928
929 if (present(atoms_moved)) atoms_moved = this%thermostat == thermo_nh
930
931 if (this%thermostat /= thermo_nh) then
932 ! velocity verlet
933
934 do iatom = 1, ions%natoms
935 if (ions%fixed(iatom)) cycle
936
937 ions%vel(:, iatom) = ions%vel(:, iatom) &
938 + this%dt/ions%mass(iatom) * m_half * (this%oldforce(:, iatom) + &
939 ions%tot_force(:, iatom))
940
941 end do
942
943 else
944 ! the nose-hoover integration
945 do iatom = 1, ions%natoms
946 ions%vel(:, iatom) = ions%vel(:, iatom) + this%dt*ions%tot_force(:, iatom) / ions%mass(iatom)
947 ions%pos(:, iatom) = ions%pos(:, iatom) + m_half*this%dt*ions%vel(:, iatom)
948 end do
949
950 call nh_chain(this, ions)
951
952 end if
953
954 if (this%thermostat == thermo_scal) then
955 scal = sqrt(this%current_temperature/ion_dynamics_temperature(ions))
956
957 ions%vel = scal*ions%vel
958 end if
959
960 if (this%cell_relax()) then
961 this%cell_vel = this%cell_vel + this%dt * m_half * (this%old_cell_force + this%cell_force)
962 end if
963
966
967
968 ! ---------------------------------------------------------
970 subroutine ion_dynamics_verlet_step1(ions, q, v, fold, dt)
971 type(ions_t), intent(in) :: ions
972 real(real64), intent(inout) :: q(:, :)
973 real(real64), intent(inout) :: v(:, :)
974 real(real64), intent(in) :: fold(:, :)
975 real(real64), intent(in) :: dt
976
977 integer :: iatom
978
980
981 ! First transform momenta to velocities
982 do iatom = 1, ions%natoms
983 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
984 end do
985
986 ! integrate using verlet
987 do iatom = 1, ions%natoms
988 if (ions%fixed(iatom)) cycle
989 q(iatom, 1:ions%space%dim) = q(iatom, 1:ions%space%dim) + dt * v(iatom, 1:ions%space%dim) + &
990 m_half*dt**2 / ions%mass(iatom) * fold(iatom, 1:ions%space%dim)
991 end do
992
993 ! And back to momenta.
994 do iatom = 1, ions%natoms
995 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
996 end do
997
999 end subroutine ion_dynamics_verlet_step1
1000
1001
1002
1003 ! ---------------------------------------------------------
1005 subroutine ion_dynamics_verlet_step2(ions, v, fold, fnew, dt)
1006 type(ions_t), intent(in) :: ions
1007 real(real64), intent(inout) :: v(:, :)
1008 real(real64), intent(in) :: fold(:, :)
1009 real(real64), intent(in) :: fnew(:, :)
1010 real(real64), intent(in) :: dt
1011
1012 integer :: iatom
1013
1015
1016 ! First transform momenta to velocities
1017 do iatom = 1, ions%natoms
1018 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
1019 end do
1020
1021 ! velocity verlet
1022 do iatom = 1, ions%natoms
1023 if (ions%fixed(iatom)) cycle
1024 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) &
1025 + dt / ions%mass(iatom) * m_half * (fold(iatom, 1:ions%space%dim) + fnew(iatom, 1:ions%space%dim))
1026 end do
1027
1028 ! And back to momenta.
1029 do iatom = 1, ions%natoms
1030 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
1031 end do
1032
1034 end subroutine ion_dynamics_verlet_step2
1035
1036
1037 ! ---------------------------------------------------------
1038 subroutine ion_dynamics_save_state(this, ions, state)
1039 type(ion_dynamics_t), intent(in) :: this
1040 type(ions_t), intent(in) :: ions
1041 type(ion_state_t), intent(out) :: state
1042
1043 if (.not. this%ions_move()) return
1044
1045 push_sub(ion_dynamics_save_state)
1046
1047 safe_allocate(state%pos(1:ions%space%dim, 1:ions%natoms))
1048 safe_allocate(state%vel(1:ions%space%dim, 1:ions%natoms))
1049
1050 state%pos = ions%pos
1051 state%vel = ions%vel
1052
1053 if (this%thermostat == thermo_nh) then
1054 safe_allocate(state%old_pos(1:ions%space%dim, 1:ions%natoms))
1055 state%old_pos(1:ions%space%dim, 1:ions%natoms) = this%old_pos(1:ions%space%dim, 1:ions%natoms)
1056 state%nh(1:2)%pos = this%nh(1:2)%pos
1057 state%nh(1:2)%vel = this%nh(1:2)%vel
1058 end if
1059
1061 end subroutine ion_dynamics_save_state
1062
1064 ! ---------------------------------------------------------
1065 subroutine ion_dynamics_restore_state(this, ions, state)
1066 type(ion_dynamics_t), intent(inout) :: this
1067 type(ions_t), intent(inout) :: ions
1068 type(ion_state_t), intent(inout) :: state
1069
1070 assert(.not. this%cell_relax())
1071 if (.not. this%ions_move()) return
1072
1074
1075 ions%pos = state%pos
1076 ions%vel = state%vel
1077
1078 if (this%thermostat == thermo_nh) then
1079 this%old_pos(1:ions%space%dim, 1:ions%natoms) = state%old_pos(1:ions%space%dim, 1:ions%natoms)
1080 this%nh(1:2)%pos = state%nh(1:2)%pos
1081 this%nh(1:2)%vel = state%nh(1:2)%vel
1082 safe_deallocate_a(state%old_pos)
1083 end if
1084
1085 safe_deallocate_a(state%pos)
1086 safe_deallocate_a(state%vel)
1087
1089 end subroutine ion_dynamics_restore_state
1090
1091
1092 ! ---------------------------------------------------------
1093 logical pure function ion_dynamics_ions_move(this) result(ions_move)
1094 class(ion_dynamics_t), intent(in) :: this
1095
1096 ions_move = this%move_ions
1097
1099
1100
1101 ! ---------------------------------------------------------
1103 logical pure function ion_dynamics_drive_ions(this) result(drive_ions)
1104 type(ion_dynamics_t), intent(in) :: this
1105
1106 drive_ions = this%drive_ions
1107
1108 end function ion_dynamics_drive_ions
1109
1110 ! ---------------------------------------------------------
1112 logical pure function ion_dynamics_cell_relax(this) result(cell_dynamics)
1113 class(ion_dynamics_t), intent(in) :: this
1114
1115 cell_dynamics = this%relax_cell
1116
1117 end function ion_dynamics_cell_relax
1118
1119 ! ---------------------------------------------------------
1121 logical pure function ion_dynamics_is_active(this) result(is_active)
1122 class(ion_dynamics_t), intent(in) :: this
1123
1124 is_active = this%relax_cell .or. this%move_ions
1125
1126 end function ion_dynamics_is_active
1127
1128
1129 ! ---------------------------------------------------------
1131 real(real64) function ion_dynamics_temperature(ions) result(temperature)
1132 type(ions_t), intent(in) :: ions
1133
1134 integer :: iatom
1135 real(real64) :: kinetic_energy
1136
1137 kinetic_energy = m_zero
1138 do iatom = 1, ions%natoms
1139 kinetic_energy = kinetic_energy + &
1140 m_half * ions%mass(iatom) * sum(ions%vel(:, iatom)**2)
1141 end do
1142 temperature = m_two/m_three*kinetic_energy/ions%natoms
1143
1144 end function ion_dynamics_temperature
1145
1146
1147 ! ---------------------------------------------------------
1149 logical function ion_dynamics_freeze(this) result(freeze)
1150 type(ion_dynamics_t), intent(inout) :: this
1151 if (this%move_ions) then
1152 this%move_ions = .false.
1153 freeze = .true.
1154 else
1155 freeze = .false.
1156 end if
1157 end function ion_dynamics_freeze
1159
1160 ! ---------------------------------------------------------
1162 subroutine ion_dynamics_unfreeze(this)
1163 type(ion_dynamics_t), intent(inout) :: this
1164 this%move_ions = .true.
1165 end subroutine ion_dynamics_unfreeze
1166
1167 ! ---------------------------------------------------------
1168 subroutine ion_dynamics_dump(this, restart, ierr)
1169 type(ion_dynamics_t), intent(in) :: this
1170 type(restart_t), intent(in) :: restart
1171 integer, intent(out) :: ierr
1172
1173 push_sub(ion_dynamics_dump)
1174
1175 if (allocated(this%oldforce)) then
1176 call drestart_write_binary(restart, "ion_dynamics_oldforce", size(this%oldforce), &
1177 this%oldforce, ierr)
1178 end if
1179
1180 if( allocated(this%old_cell_force)) then
1181 call drestart_write_binary(restart, "ion_dynamics_old_cell_force", size(this%old_cell_force), &
1182 this%old_cell_force, ierr)
1183 end if
1184
1185 pop_sub(ion_dynamics_dump)
1186 end subroutine ion_dynamics_dump
1187
1188 ! ---------------------------------------------------------
1189 subroutine ion_dynamics_load(this, restart, ierr)
1190 type(ion_dynamics_t), intent(inout) :: this
1191 type(restart_t), intent(in) :: restart
1192 integer, intent(out) :: ierr
1193
1194 push_sub(ion_dynamics_load)
1195
1196 if (allocated(this%oldforce)) then
1197 call drestart_read_binary(restart, "ion_dynamics_oldforce", size(this%oldforce), &
1198 this%oldforce, ierr)
1199 end if
1200
1201 if( allocated(this%old_cell_force)) then
1202 call drestart_read_binary(restart, "ion_dynamics_old_cell_force", size(this%old_cell_force), &
1203 this%old_cell_force, ierr)
1204 end if
1206 pop_sub(ion_dynamics_load)
1207 end subroutine ion_dynamics_load
1208
1209 !----------------------------------------------------------
1211 subroutine ion_dynamics_update_stress(this, space, stress, rlattice, rcell_volume)
1212 class(ion_dynamics_t), intent(inout) :: this
1213 class(space_t), intent(in) :: space
1214 real(real64), intent(in) :: stress(3,3)
1215 real(real64), intent(in) :: rlattice(:,:)
1216 real(real64), intent(in) :: rcell_volume
1217
1218 integer :: idir, jdir, comp
1219 real(real64) :: inv_latt(space%periodic_dim, space%periodic_dim), tmp_stress(space%periodic_dim, space%periodic_dim)
1220 real(real64) :: cell_force(space%periodic_dim, space%periodic_dim)
1221
1223
1224 ! Get the inverse lattice
1225 inv_latt = rlattice(1:space%periodic_dim, 1:space%periodic_dim)
1226 call lalg_inverse(space%periodic_dim, inv_latt, 'dir')
1227
1228 tmp_stress = -stress(1:space%periodic_dim, 1:space%periodic_dim)
1229 do idir = 1, space%periodic_dim
1230 tmp_stress(idir, idir) = tmp_stress(idir, idir) - this%pressure/space%periodic_dim
1231 end do
1232 cell_force = matmul(tmp_stress, transpose(inv_latt)) * rcell_volume
1233
1234 comp = 1
1235 do idir = 1, space%periodic_dim
1236 do jdir = 1, space%periodic_dim
1237 this%cell_force(comp) = cell_force(idir, jdir)
1238 comp = comp + 1
1239 end do
1240 end do
1241
1242 if (debug%info) then
1243 write(message(1),'(a,3a,a)') ' Stress tensor [', trim(units_abbrev(units_out%energy/units_out%length**space%dim)), ']'
1244 do idir = 1, space%periodic_dim
1245 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%energy/units_out%length**space%dim, stress(jdir, idir)), &
1246 jdir = 1, space%periodic_dim)
1247 end do
1248 call messages_info(1+space%periodic_dim, namespace=global_namespace)
1249 end if
1250
1252 end subroutine ion_dynamics_update_stress
1253
1254 !----------------------------------------------------------
1255 subroutine ion_dynamics_box_update(namespace, gr, space, new_latt)
1256 type(namespace_t), intent(in) :: namespace
1257 type(grid_t), intent(inout) :: gr
1258 class(space_t), intent(in) :: space
1259 type(lattice_vectors_t), intent(in) :: new_latt
1260
1261 integer :: idir
1262 real(real64) :: length(1:space%dim)
1263
1264 push_sub(ion_dynamics_box_update)
1265
1266 ! Regenerate the box
1267 select type(box => gr%box)
1268 type is (box_parallelepiped_t)
1269 do idir = 1, space%dim
1270 length(idir) = norm2(new_latt%rlattice(1:space%dim, idir))
1271 end do
1272 call box%regenerate(space%dim, new_latt%rlattice, length, namespace)
1273 class default
1274 call messages_not_implemented("Grid regeneration for non-parallelepiped boxes", namespace=namespace)
1275 end select
1276
1278 end subroutine ion_dynamics_box_update
1279
1280 !----------------------------------------------------------
1281 subroutine electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
1282 type(namespace_t), intent(in) :: namespace
1283 type(grid_t), intent(inout) :: gr
1284 class(space_t), intent(in) :: space
1285 type(poisson_t), intent(inout) :: psolver
1286 type(kpoints_t), intent(inout) :: kpoints
1287 type(multicomm_t), intent(in) :: mc
1288 real(real64), intent(in) :: qtot
1289 type(lattice_vectors_t), intent(in) :: new_latt
1290
1292
1293 call grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
1294
1295 !Initialize Poisson solvers
1296 call poisson_end(psolver)
1297 call poisson_init(psolver, namespace, space, gr%der, mc, gr%stencil, qtot, verbose=.false.)
1298
1299 call kpoints_lattice_vectors_update(kpoints, new_latt)
1300
1302
1305
1306end module ion_dynamics_oct_m
1307
1308!! Local Variables:
1309!! mode: f90
1310!! coding: utf-8
1311!! End:
Functions to generate random numbers.
Definition: loct_math.F90:329
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
double exp(double __x) __attribute__((__nothrow__
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public ion_dynamics_verlet_step2(ions, v, fold, fnew, dt)
A bare verlet integrator.
subroutine ion_dynamics_update_temperature(this, time, namespace)
Update the temperature of the ions in case of a thermostat.
subroutine, public ion_dynamics_dump(this, restart, ierr)
subroutine, public ion_dynamics_verlet_step1(ions, q, v, fold, dt)
A bare verlet integrator.
logical pure function ion_dynamics_cell_relax(this)
Is the cell dynamics activated or not.
subroutine nh_chain(this, ions)
subroutine, public ion_dynamics_restore_state(this, ions, state)
subroutine, public ion_dynamics_propagate(this, ions, time, dt, namespace)
Interface for the ion/cell dynamics.
subroutine ion_dynamics_propagate_cell(this, ions, time, dt, namespace)
Time-evolution of the lattice vectors.
subroutine, public ion_dynamics_save_state(this, ions, state)
subroutine, public ion_dynamics_unfreeze(this)
Unfreezes the ionic movement.
subroutine ion_dynamics_propagate_ions(this, ions, time, dt)
Time evolution of the ions.
subroutine ion_dynamics_update_stress(this, space, stress, rlattice, rcell_volume)
Updates the stress tensor for the ion dynamics.
subroutine, public ion_dynamics_propagate_vel(this, ions, atoms_moved)
subroutine, public ion_dynamics_load(this, restart, ierr)
subroutine, public ion_dynamics_init(this, namespace, ions, symmetrize, symm)
integer, parameter thermo_scal
subroutine, public electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
subroutine, public ion_dynamics_end(this)
subroutine ion_dynamics_propagate_driven_ions(this, ions, time, dt)
Move ions following a driven motion.
integer, parameter thermo_nh
logical pure function ion_dynamics_ions_move(this)
logical function, public ion_dynamics_freeze(this)
Freezes the ionic movement.
subroutine, public ion_dynamics_box_update(namespace, gr, space, new_latt)
logical pure function ion_dynamics_is_active(this)
Is the cell dynamics activated or not.
logical pure function, public ion_dynamics_drive_ions(this)
Is the ion dynamics activated or not.
real(real64) function, public ion_dynamics_temperature(ions)
This function returns the ionic temperature in energy units.
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1085
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
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 handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
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)
subroutine, public tdf_end(f)
Definition: tdfunction.F90:980
subroutine, public tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
Definition: tdfunction.F90:218
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_system_t), public units_out
type(unit_t), public unit_kelvin
For converting energies into temperatures.
type(unit_system_t), public units_inp
the units systems for reading and writing
int true(void)