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 if (ndisp > ions%natoms) call messages_input_error(namespace, 'IonsTimeDependentDisplacements')
278
279 do i = 1, ndisp
280 call parse_block_integer(blk, i-1, 0, iatom)
281 this%td_displacements(iatom)%move = .true.
282 if (iatom < 1 .and. iatom > ions%natoms) then
283 call messages_input_error(namespace, 'IonsTimeDependentDisplacements')
284 end if
285
286 call parse_block_string(blk, i-1, 1, expression)
287 call tdf_read(this%td_displacements(iatom)%fx, namespace, trim(expression), ierr)
288 if (ierr /= 0) then
289 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
290 call messages_warning(1, namespace=namespace)
291 end if
292
293
294 call parse_block_string(blk, i-1, 2, expression)
295 call tdf_read(this%td_displacements(iatom)%fy, namespace, trim(expression), ierr)
296 if (ierr /= 0) then
297 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
298 call messages_warning(1, namespace=namespace)
299 end if
300
301 call parse_block_string(blk, i-1, 3, expression)
302 call tdf_read(this%td_displacements(iatom)%fz, namespace, trim(expression), ierr)
303 if (ierr /= 0) then
304 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
305 call messages_warning(1, namespace=namespace)
306 end if
307
308 end do
309
310 safe_allocate(this%ions_t0)
311 this%ions_t0 = ions
312
313 end if
314
315
316
317
318 !%Variable Thermostat
319 !%Type integer
320 !%Default none
321 !%Section Time-Dependent::Propagation
322 !%Description
323 !% This variable selects the type of thermostat applied to
324 !% control the ionic temperature.
325 !%Option none 0
326 !% No thermostat is applied. This is the default.
327 !%Option velocity_scaling 1
328 !% Velocities are scaled to control the temperature.
329 !%Option nose_hoover 2
330 !% Nose-Hoover thermostat.
331 !%End
332
333 call parse_variable(namespace, 'Thermostat', thermo_none, this%thermostat)
334 if (.not. varinfo_valid_option('Thermostat', this%thermostat)) call messages_input_error(namespace, 'Thermostat')
335 call messages_print_var_option('Thermostat', this%thermostat, namespace=namespace)
336
337 if (this%thermostat /= thermo_none) then
338
339 have_velocities = .true.
340
341 if (this%drive_ions) then
342 call messages_write('You cannot use a Thermostat and IonsConstantVelocity or IonsTimeDependentDisplacements')
343 call messages_write('at the same time.')
344 call messages_fatal(namespace=namespace)
345 end if
346
347 call messages_experimental('Thermostat', namespace=namespace)
348
349 !%Variable TemperatureFunction
350 !%Type integer
351 !%Default "temperature"
352 !%Section Time-Dependent::Propagation
353 !%Description
354 !% If a thermostat is used, this variable indicates the name of the
355 !% function in the <tt>TDFunctions</tt> block that will be used to control the
356 !% temperature. The values of the temperature are given in
357 !% degrees Kelvin.
358 !%End
359 call parse_variable(namespace, 'TemperatureFunction', 'temperature', temp_function_name)
360
361 call tdf_read(this%temperature_function, namespace, temp_function_name, ierr)
362
363 if (ierr /= 0) then
364 message(1) = "You have enabled a thermostat but Octopus could not find"
365 message(2) = "the '"//trim(temp_function_name)//"' function in the TDFunctions block."
366 call messages_fatal(2, namespace=namespace)
367 end if
368
369 if (this%thermostat == thermo_nh) then
370 !%Variable ThermostatMass
371 !%Type float
372 !%Default 1.0
373 !%Section Time-Dependent::Propagation
374 !%Description
375 !% This variable sets the fictitious mass for the Nose-Hoover
376 !% thermostat.
377 !%End
378 call messages_obsolete_variable(namespace, 'NHMass', 'ThermostatMass')
379
380 call parse_variable(namespace, 'ThermostatMass', m_one, this%nh(1)%mass)
381 this%nh(2)%mass = this%nh(1)%mass
382
383 this%nh(1:2)%pos = m_zero
384 this%nh(1:2)%vel = m_zero
385
386 safe_allocate(this%old_pos(1:ions%space%dim, 1:ions%natoms))
387
388 this%old_pos = ions%pos
389 end if
390
391 end if
392
393 !now initialize velocities
394
395 !%Variable RandomVelocityTemp
396 !%Type float
397 !%Default 0.0
398 !%Section System::Velocities
399 !%Description
400 !% If this variable is present, <tt>Octopus</tt> will assign random
401 !% velocities to the atoms following a Boltzmann distribution with
402 !% temperature given by <tt>RandomVelocityTemp</tt> (in degrees Kelvin).
403 !% The seed for the random number generator can be modified by setting
404 !% <tt>GSL_RNG_SEED</tt> environment variable.
405 !%End
406
407 ! we now load the velocities, either from the temperature, from the input, or from a file
408 if (parse_is_defined(namespace, 'RandomVelocityTemp')) then
409
410 have_velocities = .true.
411
412 if (mpi_grp_is_root(mpi_world)) then
413 call loct_ran_init(random_gen_pointer)
414 call parse_variable(namespace, 'RandomVelocityTemp', m_zero, temperature, unit = unit_kelvin)
415 end if
416
417 do i = 1, ions%natoms
418 !generate the velocities in the root node
419 if (mpi_grp_is_root(mpi_world)) then
420 sigma = sqrt(temperature / ions%mass(i))
421 do j = 1, 3
422 ions%vel(j, i) = loct_ran_gaussian(random_gen_pointer, sigma)
423 end do
424 end if
425 !and send them to the others
426 call mpi_world%bcast(ions%vel(:, i), ions%space%dim, mpi_double_precision, 0)
427 end do
428
429 if (mpi_grp_is_root(mpi_world)) then
430 call loct_ran_end(random_gen_pointer)
431 end if
432
433 call ions%update_kinetic_energy()
434 kin1 = ions%kinetic_energy
435
436 xx = ions%center_of_mass_vel()
437 do i = 1, ions%natoms
438 ions%vel(:, i) = ions%vel(:, i) - xx
439 end do
440
441 call ions%update_kinetic_energy()
442 kin2 = ions%kinetic_energy
443
444 if (kin2 > m_epsilon) then
445 do i = 1, ions%natoms
446 ions%vel(:, i) = sqrt(kin1/kin2)*ions%vel(:, i)
447 end do
448 end if
449
450 call ions%update_kinetic_energy()
451
452 write(message(1),'(a,f10.4,1x,a)') 'Info: Initial velocities randomly distributed with T =', &
454 write(message(2),'(2x,a,f8.4,1x,a)') '<K> =', &
455 units_from_atomic(units_out%energy, ions%kinetic_energy/ions%natoms), &
456 units_abbrev(units_out%energy)
457 write(message(3),'(2x,a,f8.4,1x,a)') '3/2 k_B T =', &
458 units_from_atomic(units_out%energy, (m_three/m_two)*temperature), &
459 units_abbrev(units_out%energy)
460 call messages_info(3, namespace=namespace)
461
462 else
463 !%Variable XYZVelocities
464 !%Type string
465 !%Section System::Velocities
466 !%Description
467 !% <tt>Octopus</tt> will try to read the starting velocities of the atoms from the XYZ file
468 !% specified by the variable <tt>XYZVelocities</tt>.
469 !% Note that you do not need to specify initial velocities if you are not going
470 !% to perform ion dynamics; if you are going to allow the ions to move but the velocities
471 !% are not specified, they are considered to be null.
472 !% Note: It is important for the velocities to maintain the ordering
473 !% in which the atoms were defined in the coordinates specifications.
474 !%End
475
476 !%Variable XSFVelocities
477 !%Type string
478 !%Section System::Velocities
479 !%Description
480 !% Like <tt>XYZVelocities</tt> but in XCrySDen format, as in <tt>XSFCoordinates</tt>.
481 !%End
482
483 !%Variable PDBVelocities
484 !%Type string
485 !%Section System::Velocities
486 !%Description
487 !% Like <tt>XYZVelocities</tt> but in PDB format, as in <tt>PDBCoordinates</tt>.
488 !%End
489
490 !%Variable Velocities
491 !%Type block
492 !%Section System::Velocities
493 !%Description
494 !% If <tt>XYZVelocities</tt>, <tt>PDBVelocities</tt>, and <tt>XSFVelocities</tt>
495 !% are not present, <tt>Octopus</tt> will try to fetch the initial
496 !% atomic velocities from this block. If this block is not present, <tt>Octopus</tt>
497 !% will set the initial velocities to zero. The format of this block can be
498 !% illustrated by this example:
499 !%
500 !% <tt>%Velocities
501 !% <br>&nbsp;&nbsp;'C' | -1.7 | 0.0 | 0.0
502 !% <br>&nbsp;&nbsp;'O' | &nbsp;1.7 | 0.0 | 0.0
503 !% <br>%</tt>
504 !%
505 !% It describes one carbon and one oxygen moving at the relative
506 !% velocity of 3.4 velocity units.
507 !%
508 !% Note: It is important for the velocities to maintain the ordering
509 !% in which the atoms were defined in the coordinates specifications.
510 !%End
511
512 call read_coords_init(xyz)
513 call read_coords_read('Velocities', xyz, ions%space, namespace)
514 if (xyz%source /= read_coords_err) then
515
516 have_velocities = .true.
517
518 if (ions%natoms /= xyz%n) then
519 write(message(1), '(a,i4,a,i4)') 'I need exactly ', ions%natoms, ' velocities, but I found ', xyz%n
520 call messages_fatal(1, namespace=namespace)
521 end if
522
523 ! copy information and adjust units
524 do i = 1, ions%natoms
525 ions%vel(:, i) = units_to_atomic(units_inp%velocity/units_inp%length, xyz%atom(i)%x(1:ions%space%dim))
526 end do
527
528 call read_coords_end(xyz)
529
530 else
531 ions%vel = m_zero
532 end if
533 end if
534
535 call ions%update_kinetic_energy()
536
537 !%Variable MoveIons
538 !%Type logical
539 !%Section Time-Dependent::Propagation
540 !%Description
541 !% This variable controls whether atoms are moved during a time
542 !% propagation run. The default is yes when the ion velocity is
543 !% set explicitly or implicitly, otherwise is no.
544 !%End
545 call parse_variable(namespace, 'MoveIons', have_velocities, this%move_ions)
546 call messages_print_var_value('MoveIons', this%move_ions, namespace=namespace)
547
548 if (this%move_ions .and. ions%space%periodic_dim == 1) then
549 call messages_input_error(namespace, 'MoveIons', &
550 'Moving ions for a 1D periodic system is not allowed, as forces are incorrect.')
551 end if
552
553 if (this%ions_move()) then
554 safe_allocate(this%oldforce(1:ions%space%dim, 1:ions%natoms))
555 end if
556
557 if (ions%space%is_periodic()) then
558 !%Variable CellDynamics
559 !%Type logical
560 !%Section Time-Dependent::Propagation
561 !%Description
562 !% This variable controls whether the cell relaxation is done during a time
563 !% propagation run. The default is no.
564 !% This is done based on the Parrinello-Rahman equation of motion of the cell,
565 !% see Parrinello and Rahman, J. Appl. Pys. 52, 7182 (1981).
566 !% Only for periodic systems.
567 !%End
568 call parse_variable(namespace, 'CellDynamics', .false., this%relax_cell)
569 call messages_print_var_value('CellDynamics', this%relax_cell, namespace=namespace)
570
571 if (this%cell_relax()) then
572 if (accel_is_enabled()) then
573 message(1) = "Cell dynamics not supported on GPUs."
574 call messages_fatal(1, namespace=namespace)
575 end if
576
577 periodic_dim = ions%space%periodic_dim
578 ncomp = periodic_dim * periodic_dim
579 safe_allocate(this%cell_force(1:ncomp))
580 this%cell_force = m_zero
581 safe_allocate(this%old_cell_force(1:ncomp))
582 this%old_cell_force = m_zero
583 safe_allocate(this%cell_vel(1:ncomp))
584 this%cell_vel = m_zero
585 ! We start from identity for the strain
586 safe_allocate(this%strain(1:ncomp))
587 this%strain = m_zero
588 ncomp = 1
589 do i = 1, periodic_dim
590 do j = i, periodic_dim
591 if(i == j) this%strain(ncomp) = m_one
592 ncomp = ncomp + 1
593 end do
594 end do
595
596 ! As we work with "reduced" quantities, the initial lattice vectors are used as a reference
597 ! and the strain is propagated from this initial reference
598 safe_allocate(this%initial_rlattice(1:periodic_dim, 1:periodic_dim))
599 this%initial_rlattice(1:periodic_dim, 1:periodic_dim) = ions%latt%rlattice(1:periodic_dim, 1:periodic_dim)
600 end if
601
602 !%Variable HydrostaticPressure
603 !%Type float
604 !%Default 0.0
605 !%Section Time-Dependent::Propagation
606 !%Description
607 !% Geometry optimization and molecular dynamics can be performed in presence of an external
608 !% hydrostatic pressure. This variable allows to set this value.
609 !% Only for periodic systems.
610 !%End
611 call parse_variable(namespace, 'HydrostaticPressure', m_zero, this%pressure)
612 else
613 this%relax_cell = .false.
614 end if
615
616
617 pop_sub(ion_dynamics_init)
618 end subroutine ion_dynamics_init
619
620
621 ! ---------------------------------------------------------
622 subroutine ion_dynamics_end(this)
623 type(ion_dynamics_t), intent(inout) :: this
624
625 push_sub(ion_dynamics_end)
626 safe_deallocate_a(this%oldforce)
627
628 if (this%thermostat /= thermo_none) then
629 call tdf_end(this%temperature_function)
630 end if
631
632 if (this%drive_ions .and. allocated(this%td_displacements)) then
633 if (any(this%td_displacements(1:this%ions_t0%natoms)%move)) then
634 ! ions end cannot be called here, otherwise the species are destroyed twice
635 ! call ions_end(this%ions_t0)
636 !AFE_DEALLOCATE_P(this%ions_t0)
637 end if
638 safe_deallocate_a(this%td_displacements)
639 end if
640
641 safe_deallocate_a(this%cell_force)
642 safe_deallocate_a(this%old_cell_force)
643 safe_deallocate_a(this%cell_vel)
644 safe_deallocate_a(this%initial_rlattice)
645
646 pop_sub(ion_dynamics_end)
647 end subroutine ion_dynamics_end
648
649
650 ! ---------------------------------------------------------
652 subroutine ion_dynamics_propagate(this, ions, time, dt, namespace)
653 type(ion_dynamics_t), intent(inout) :: this
654 type(ions_t), intent(inout) :: ions
655 real(real64), intent(in) :: time
656 real(real64), intent(in) :: dt
657 type(namespace_t), intent(in) :: namespace
658
659 integer :: iatom
660
661 push_sub(ion_dynamics_propagate)
662
663 this%dt = dt
664
665 if (this%drive_ions) then
666
667 call ion_dynamics_propagate_driven_ions(this, ions, time, dt)
668
669 else
670
671 ! Update the thermostat temperature
672 call ion_dynamics_update_temperature(this, time, namespace)
673
674 ! Conversion to reduced coordinates
675 do iatom = 1, ions%natoms
676 ions%pos(:, iatom) = ions%latt%cart_to_red(ions%pos(:, iatom))
677 ions%vel(:, iatom) = ions%latt%cart_to_red(ions%vel(:, iatom))
678 ions%tot_force(:, iatom) = ions%latt%cart_to_red(ions%tot_force(:, iatom))
679 end do
680
681 ! Get the new reduced coordinates
682 if (this%ions_move()) then
683 call ion_dynamics_propagate_ions(this, ions, time, dt)
684 end if
685
686 ! Updating the lattice vectors
687 if (this%cell_relax()) then
688 call ion_dynamics_propagate_cell(this, ions, time, dt, namespace)
689 end if
690
691 ! Get the new Cartesian coordinates
692 do iatom = 1, ions%natoms
693 ions%pos(:, iatom) = ions%latt%red_to_cart(ions%pos(:, iatom))
694 ions%vel(:, iatom) = ions%latt%red_to_cart(ions%vel(:, iatom))
695 if (allocated(this%oldforce)) then
696 this%oldforce(:, iatom) = ions%latt%red_to_cart(this%oldforce(:, iatom))
697 end if
698 ions%tot_force(:, iatom) = ions%latt%red_to_cart(ions%tot_force(:, iatom))
699 end do
700
701 end if
702
703 call ions%fold_atoms_into_cell()
704
706 end subroutine ion_dynamics_propagate
707
708
709 ! ---------------------------------------------------------
711 subroutine ion_dynamics_update_temperature(this, time, namespace)
712 type(ion_dynamics_t), intent(inout) :: this
713 real(real64), intent(in) :: time
714 type(namespace_t), intent(in) :: namespace
717
718 ! get the temperature from the tdfunction for the current time
719 if (this%thermostat /= thermo_none) then
720 this%current_temperature = units_to_atomic(unit_kelvin, tdf(this%temperature_function, time))
721
722 if (this%current_temperature < m_zero) then
723 write(message(1), '(a, f10.3, 3a, f10.3, 3a)') &
724 "Negative temperature (", &
725 units_from_atomic(unit_kelvin, this%current_temperature), " ", units_abbrev(unit_kelvin), &
726 ") at time ", &
727 units_from_atomic(units_out%time, time), " ", trim(units_abbrev(units_out%time)), "."
728 call messages_fatal(1, namespace=namespace)
729 end if
730 else
731 this%current_temperature = m_zero
732 end if
733
736
740 subroutine ion_dynamics_propagate_driven_ions(this, ions, time, dt)
741 type(ion_dynamics_t), intent(inout) :: this
742 type(ions_t), intent(inout) :: ions
743 real(real64), intent(in) :: time
744 real(real64), intent(in) :: dt
746 integer :: iatom
747 real(real64) :: dr(3)
748
750
751 assert(this%drive_ions)
752
753 do iatom = 1, ions%natoms
754 if (ions%fixed(iatom)) cycle
755
756 if (this%constant_velocity) then
757 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom)
758
759 else if (allocated(this%td_displacements)) then
760
761 if (this%td_displacements(iatom)%move) then
762 dr(1:3)=(/ real(tdf(this%td_displacements(iatom)%fx, time), real64), &
763 real(tdf(this%td_displacements(iatom)%fy, time), real64), &
764 real(tdf(this%td_displacements(iatom)%fz, time), real64) /)
765
766 ions%pos(:, iatom) = this%ions_t0%pos(:, iatom) + dr(1:ions%space%dim)
767 end if
768
769 end if
770 end do
771
774
775 ! ---------------------------------------------------------
780 subroutine ion_dynamics_propagate_ions(this, ions, time, dt)
781 type(ion_dynamics_t), intent(inout) :: this
782 type(ions_t), intent(inout) :: ions
783 real(real64), intent(in) :: time
784 real(real64), intent(in) :: dt
785
786 integer :: iatom
787
789
790 assert(.not. this%drive_ions)
791
792 if (this%thermostat /= thermo_nh) then
793 ! integrate using verlet
794 do iatom = 1, ions%natoms
795 if (ions%fixed(iatom)) cycle
796
797 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom) + &
798 m_half*dt**2 / ions%mass(iatom) * ions%tot_force(:, iatom)
799
800 this%oldforce(:, iatom) = ions%tot_force(:, iatom)
801 end do
802
803 else
804 ! for the Nose-Hoover thermostat we use a special integrator
805
806 ! The implementation of the Nose-Hoover thermostat is based on
807 ! Understanding Molecular Simulations by Frenkel and Smit,
808 ! Appendix E, page 540-542.
809
810 call nh_chain(this, ions)
811
812 do iatom = 1, ions%natoms
813 if (ions%fixed(iatom)) cycle
814
815 ions%pos(:, iatom) = ions%pos(:, iatom) + m_half*dt*ions%vel(:, iatom)
816 end do
817
818 end if
819
821 end subroutine ion_dynamics_propagate_ions
822
823 ! ---------------------------------------------------------
825 subroutine ion_dynamics_propagate_cell(this, ions, time, dt, namespace)
826 type(ion_dynamics_t), intent(inout) :: this
827 type(ions_t), intent(inout) :: ions
828 real(real64), intent(in) :: time
829 real(real64), intent(in) :: dt
830 type(namespace_t), intent(in) :: namespace
831
832 integer :: idir, jdir, comp
833 real(real64) :: rlattice_change(ions%space%periodic_dim*ions%space%periodic_dim)
834
836
837 rlattice_change = dt * this%cell_vel + m_half*dt**2 * this%cell_force
838
839 write(message(1),'(a,3a,a)') ' Cell force [', trim(units_abbrev(units_out%energy/units_out%length**ions%space%dim)), ']'
840 do idir = 1, ions%space%periodic_dim
841 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%energy/units_out%length**ions%space%dim, &
842 this%cell_force(jdir + (idir-1)*ions%space%periodic_dim)), &
843 jdir = 1, ions%space%periodic_dim)
844 end do
845 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
846
847 write(message(1),'(a,3a,a)') ' Cell vel [', &
848 trim(units_abbrev(units_out%energy/units_out%length**ions%space%dim*units_out%time)), ']'
849 do idir = 1, ions%space%periodic_dim
850 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%energy/units_out%length**ions%space%dim*units_out%time, &
851 this%cell_vel(jdir+ (idir-1)*ions%space%periodic_dim)), &
852 jdir = 1, ions%space%periodic_dim)
853 end do
854 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
855
856
857 comp = 1
858 do idir = 1, ions%space%periodic_dim
859 do jdir = 1, ions%space%periodic_dim
860 ions%latt%rlattice(idir, jdir) = ions%latt%rlattice(idir, jdir) + rlattice_change(comp)
861 comp = comp + 1
862 end do
863 end do
864
865 this%old_cell_force = this%cell_force
866
867 if (associated(this%symm)) then
868 call this%symm%symmetrize_lattice_vectors(ions%space%periodic_dim, &
869 this%initial_rlattice, ions%latt%rlattice, this%symmetrize)
870 end if
871 call ions%update_lattice_vectors(ions%latt, this%symmetrize)
872
874 end subroutine ion_dynamics_propagate_cell
875
876
877 ! ---------------------------------------------------------
878 subroutine nh_chain(this, ions)
879 type(ion_dynamics_t), intent(inout) :: this
880 type(ions_t), intent(inout) :: ions
881
882 real(real64) :: g1, g2, ss, uk, dt, temp
883
884 push_sub(nh_chain)
885
886 dt = this%dt
887
888 call ions%update_kinetic_energy()
889 uk = ions%kinetic_energy
890
891 temp = this%current_temperature
892
893 g2 = (this%nh(1)%mass*this%nh(1)%vel**2 - temp)/this%nh(2)%mass
894 this%nh(2)%vel = this%nh(2)%vel + g2*dt/m_four
895 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
896
897 g1 = (m_two*uk - m_three*ions%natoms*temp)/this%nh(1)%mass
898 this%nh(1)%vel = this%nh(1)%vel + g1*dt/m_four
899 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
900 this%nh(1)%pos = this%nh(1)%pos + this%nh(1)%vel*dt/m_two
901 this%nh(2)%pos = this%nh(2)%pos + this%nh(2)%vel*dt/m_two
902
903 ss = exp(-this%nh(1)%vel*dt/m_two)
904
905 ions%vel = ss*ions%vel
906
907 uk = uk*ss**2
908
909 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
910 g1 = (m_two*uk - m_three*ions%natoms*temp)/this%nh(1)%mass
911 this%nh(1)%vel = this%nh(1)%vel + g1*dt/m_four
912 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
913
914 g2 = (this%nh(1)%mass*this%nh(1)%vel**2 - temp)/this%nh(2)%mass
915 this%nh(2)%vel = this%nh(2)%vel + g2*dt/m_four
916
917 pop_sub(nh_chain)
918 end subroutine nh_chain
919
920
921 ! ---------------------------------------------------------
922 subroutine ion_dynamics_propagate_vel(this, ions, atoms_moved)
923 type(ion_dynamics_t), intent(inout) :: this
924 type(ions_t), intent(inout) :: ions
925 logical, optional, intent(out) :: atoms_moved
926
927 integer :: iatom
928 real(real64) :: scal, temp
929
930 if (.not. ion_dynamics_ions_move(this)) return
931 if (this%drive_ions) return
932
934
935 if (present(atoms_moved)) atoms_moved = this%thermostat == thermo_nh
936
937 if (this%thermostat /= thermo_nh) then
938 ! velocity verlet
939
940 do iatom = 1, ions%natoms
941 if (ions%fixed(iatom)) cycle
942
943 ions%vel(:, iatom) = ions%vel(:, iatom) &
944 + this%dt/ions%mass(iatom) * m_half * (this%oldforce(:, iatom) + &
945 ions%tot_force(:, iatom))
946
947 end do
948
949 else
950 ! the nose-hoover integration
951 do iatom = 1, ions%natoms
952 if (ions%fixed(iatom)) cycle
953
954 ions%vel(:, iatom) = ions%vel(:, iatom) + this%dt*ions%tot_force(:, iatom) / ions%mass(iatom)
955 ions%pos(:, iatom) = ions%pos(:, iatom) + m_half*this%dt*ions%vel(:, iatom)
956 end do
957
958 call nh_chain(this, ions)
959
960 end if
961
962 if (this%thermostat == thermo_scal) then
963 temp = ion_dynamics_temperature(ions)
964 if (temp > m_epsilon) then
965 scal = sqrt(this%current_temperature/temp)
966 ions%vel = scal*ions%vel
967 end if
968 end if
969
970 if (this%cell_relax()) then
971 this%cell_vel = this%cell_vel + this%dt * m_half * (this%old_cell_force + this%cell_force)
972 end if
973
975 end subroutine ion_dynamics_propagate_vel
976
977
978 ! ---------------------------------------------------------
982 subroutine ion_dynamics_verlet_step1(ions, q, v, fold, dt)
983 type(ions_t), intent(in) :: ions
984 real(real64), intent(inout) :: q(:, :)
985 real(real64), intent(inout) :: v(:, :)
986 real(real64), intent(in) :: fold(:, :)
987 real(real64), intent(in) :: dt
988
989 integer :: iatom
990
992
993 ! First transform momenta to velocities
994 do iatom = 1, ions%natoms
995 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
996 end do
997
998 ! integrate using verlet
999 do iatom = 1, ions%natoms
1000 if (ions%fixed(iatom)) cycle
1001 q(iatom, 1:ions%space%dim) = q(iatom, 1:ions%space%dim) + dt * v(iatom, 1:ions%space%dim) + &
1002 m_half*dt**2 / ions%mass(iatom) * fold(iatom, 1:ions%space%dim)
1003 end do
1004
1005 ! And back to momenta.
1006 do iatom = 1, ions%natoms
1007 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
1008 end do
1009
1011 end subroutine ion_dynamics_verlet_step1
1012
1013
1014
1015 ! ---------------------------------------------------------
1017 subroutine ion_dynamics_verlet_step2(ions, v, fold, fnew, dt)
1018 type(ions_t), intent(in) :: ions
1019 real(real64), intent(inout) :: v(:, :)
1020 real(real64), intent(in) :: fold(:, :)
1021 real(real64), intent(in) :: fnew(:, :)
1022 real(real64), intent(in) :: dt
1023
1024 integer :: iatom
1025
1027
1028 ! First transform momenta to velocities
1029 do iatom = 1, ions%natoms
1030 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
1031 end do
1032
1033 ! velocity verlet
1034 do iatom = 1, ions%natoms
1035 if (ions%fixed(iatom)) cycle
1036 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) &
1037 + dt / ions%mass(iatom) * m_half * (fold(iatom, 1:ions%space%dim) + fnew(iatom, 1:ions%space%dim))
1038 end do
1039
1040 ! And back to momenta.
1041 do iatom = 1, ions%natoms
1042 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
1043 end do
1044
1046 end subroutine ion_dynamics_verlet_step2
1047
1048
1049 ! ---------------------------------------------------------
1050 subroutine ion_dynamics_save_state(this, ions, state)
1051 type(ion_dynamics_t), intent(in) :: this
1052 type(ions_t), intent(in) :: ions
1053 type(ion_state_t), intent(out) :: state
1054
1055 if (.not. this%ions_move()) return
1056
1057 push_sub(ion_dynamics_save_state)
1058
1059 safe_allocate(state%pos(1:ions%space%dim, 1:ions%natoms))
1060 safe_allocate(state%vel(1:ions%space%dim, 1:ions%natoms))
1061
1062 state%pos = ions%pos
1063 state%vel = ions%vel
1064
1065 if (this%thermostat == thermo_nh) then
1066 safe_allocate(state%old_pos(1:ions%space%dim, 1:ions%natoms))
1067 state%old_pos(1:ions%space%dim, 1:ions%natoms) = this%old_pos(1:ions%space%dim, 1:ions%natoms)
1068 state%nh(1:2)%pos = this%nh(1:2)%pos
1069 state%nh(1:2)%vel = this%nh(1:2)%vel
1070 end if
1071
1073 end subroutine ion_dynamics_save_state
1074
1076 ! ---------------------------------------------------------
1077 subroutine ion_dynamics_restore_state(this, ions, state)
1078 type(ion_dynamics_t), intent(inout) :: this
1079 type(ions_t), intent(inout) :: ions
1080 type(ion_state_t), intent(inout) :: state
1081
1082 assert(.not. this%cell_relax())
1083 if (.not. this%ions_move()) return
1084
1086
1087 ions%pos = state%pos
1088 ions%vel = state%vel
1089
1090 if (this%thermostat == thermo_nh) then
1091 this%old_pos(1:ions%space%dim, 1:ions%natoms) = state%old_pos(1:ions%space%dim, 1:ions%natoms)
1092 this%nh(1:2)%pos = state%nh(1:2)%pos
1093 this%nh(1:2)%vel = state%nh(1:2)%vel
1094 safe_deallocate_a(state%old_pos)
1095 end if
1096
1097 safe_deallocate_a(state%pos)
1098 safe_deallocate_a(state%vel)
1099
1101 end subroutine ion_dynamics_restore_state
1102
1103
1104 ! ---------------------------------------------------------
1105 logical pure function ion_dynamics_ions_move(this) result(ions_move)
1106 class(ion_dynamics_t), intent(in) :: this
1107
1108 ions_move = this%move_ions
1109
1111
1112
1113 ! ---------------------------------------------------------
1115 logical pure function ion_dynamics_drive_ions(this) result(drive_ions)
1116 type(ion_dynamics_t), intent(in) :: this
1117
1118 drive_ions = this%drive_ions
1119
1120 end function ion_dynamics_drive_ions
1121
1122 ! ---------------------------------------------------------
1124 logical pure function ion_dynamics_cell_relax(this) result(cell_dynamics)
1125 class(ion_dynamics_t), intent(in) :: this
1126
1127 cell_dynamics = this%relax_cell
1128
1129 end function ion_dynamics_cell_relax
1130
1131 ! ---------------------------------------------------------
1133 logical pure function ion_dynamics_is_active(this) result(is_active)
1134 class(ion_dynamics_t), intent(in) :: this
1135
1136 is_active = this%relax_cell .or. this%move_ions
1137
1138 end function ion_dynamics_is_active
1139
1140
1141 ! ---------------------------------------------------------
1143 real(real64) function ion_dynamics_temperature(ions) result(temperature)
1144 type(ions_t), intent(in) :: ions
1145
1146 integer :: iatom
1147 real(real64) :: kinetic_energy
1148
1149 kinetic_energy = m_zero
1150 do iatom = 1, ions%natoms
1151 kinetic_energy = kinetic_energy + &
1152 m_half * ions%mass(iatom) * sum(ions%vel(:, iatom)**2)
1153 end do
1154 temperature = m_two/m_three*kinetic_energy/ions%natoms
1155
1156 end function ion_dynamics_temperature
1157
1158
1159 ! ---------------------------------------------------------
1161 logical function ion_dynamics_freeze(this) result(freeze)
1162 type(ion_dynamics_t), intent(inout) :: this
1163 if (this%move_ions) then
1164 this%move_ions = .false.
1165 freeze = .true.
1166 else
1167 freeze = .false.
1168 end if
1169 end function ion_dynamics_freeze
1171
1172 ! ---------------------------------------------------------
1174 subroutine ion_dynamics_unfreeze(this)
1175 type(ion_dynamics_t), intent(inout) :: this
1176 this%move_ions = .true.
1177 end subroutine ion_dynamics_unfreeze
1178
1179 ! ---------------------------------------------------------
1180 subroutine ion_dynamics_dump(this, restart, ierr)
1181 type(ion_dynamics_t), intent(in) :: this
1182 type(restart_t), intent(in) :: restart
1183 integer, intent(out) :: ierr
1184
1185 push_sub(ion_dynamics_dump)
1186
1187 if (allocated(this%oldforce)) then
1188 call drestart_write_binary(restart, "ion_dynamics_oldforce", size(this%oldforce), &
1189 this%oldforce, ierr)
1190 end if
1191
1192 if( allocated(this%old_cell_force)) then
1193 call drestart_write_binary(restart, "ion_dynamics_old_cell_force", size(this%old_cell_force), &
1194 this%old_cell_force, ierr)
1195 end if
1196
1197 pop_sub(ion_dynamics_dump)
1198 end subroutine ion_dynamics_dump
1199
1200 ! ---------------------------------------------------------
1201 subroutine ion_dynamics_load(this, restart, ierr)
1202 type(ion_dynamics_t), intent(inout) :: this
1203 type(restart_t), intent(in) :: restart
1204 integer, intent(out) :: ierr
1205
1206 push_sub(ion_dynamics_load)
1207
1208 if (allocated(this%oldforce)) then
1209 call drestart_read_binary(restart, "ion_dynamics_oldforce", size(this%oldforce), &
1210 this%oldforce, ierr)
1211 end if
1212
1213 if( allocated(this%old_cell_force)) then
1214 call drestart_read_binary(restart, "ion_dynamics_old_cell_force", size(this%old_cell_force), &
1215 this%old_cell_force, ierr)
1216 end if
1218 pop_sub(ion_dynamics_load)
1219 end subroutine ion_dynamics_load
1220
1221 !----------------------------------------------------------
1223 subroutine ion_dynamics_update_stress(this, space, stress, rlattice, rcell_volume)
1224 class(ion_dynamics_t), intent(inout) :: this
1225 class(space_t), intent(in) :: space
1226 real(real64), intent(in) :: stress(3,3)
1227 real(real64), intent(in) :: rlattice(:,:)
1228 real(real64), intent(in) :: rcell_volume
1229
1230 integer :: idir, jdir, comp
1231 real(real64) :: inv_latt(space%periodic_dim, space%periodic_dim), tmp_stress(space%periodic_dim, space%periodic_dim)
1232 real(real64) :: cell_force(space%periodic_dim, space%periodic_dim)
1233
1235
1236 ! Get the inverse lattice
1237 inv_latt = rlattice(1:space%periodic_dim, 1:space%periodic_dim)
1238 call lalg_inverse(space%periodic_dim, inv_latt, 'dir')
1239
1240 tmp_stress = -stress(1:space%periodic_dim, 1:space%periodic_dim)
1241 do idir = 1, space%periodic_dim
1242 tmp_stress(idir, idir) = tmp_stress(idir, idir) - this%pressure/space%periodic_dim
1243 end do
1244 cell_force = matmul(tmp_stress, transpose(inv_latt)) * rcell_volume
1245
1246 comp = 1
1247 do idir = 1, space%periodic_dim
1248 do jdir = 1, space%periodic_dim
1249 this%cell_force(comp) = cell_force(idir, jdir)
1250 comp = comp + 1
1251 end do
1252 end do
1253
1254 if (debug%info) then
1255 write(message(1),'(a,3a,a)') ' Stress tensor [', trim(units_abbrev(units_out%energy/units_out%length**space%dim)), ']'
1256 do idir = 1, space%periodic_dim
1257 write(message(1+idir),'(9e18.10)') (units_from_atomic(units_out%energy/units_out%length**space%dim, stress(jdir, idir)), &
1258 jdir = 1, space%periodic_dim)
1259 end do
1260 call messages_info(1+space%periodic_dim, namespace=global_namespace)
1261 end if
1262
1264 end subroutine ion_dynamics_update_stress
1265
1266 !----------------------------------------------------------
1267 subroutine ion_dynamics_box_update(namespace, gr, space, new_latt)
1268 type(namespace_t), intent(in) :: namespace
1269 type(grid_t), intent(inout) :: gr
1270 class(space_t), intent(in) :: space
1271 type(lattice_vectors_t), intent(in) :: new_latt
1272
1273 integer :: idir
1274 real(real64) :: length(1:space%dim)
1275
1276 push_sub(ion_dynamics_box_update)
1277
1278 ! Regenerate the box
1279 select type(box => gr%box)
1280 type is (box_parallelepiped_t)
1281 do idir = 1, space%dim
1282 length(idir) = norm2(new_latt%rlattice(1:space%dim, idir))
1283 end do
1284 call box%regenerate(space%dim, new_latt%rlattice, length, namespace)
1285 class default
1286 call messages_not_implemented("Grid regeneration for non-parallelepiped boxes", namespace=namespace)
1287 end select
1288
1290 end subroutine ion_dynamics_box_update
1291
1292 !----------------------------------------------------------
1293 subroutine electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
1294 type(namespace_t), intent(in) :: namespace
1295 type(grid_t), intent(inout) :: gr
1296 class(space_t), intent(in) :: space
1297 type(poisson_t), intent(inout) :: psolver
1298 type(kpoints_t), intent(inout) :: kpoints
1299 type(multicomm_t), intent(in) :: mc
1300 real(real64), intent(in) :: qtot
1301 type(lattice_vectors_t), intent(in) :: new_latt
1302
1304
1305 call grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
1306
1307 !Initialize Poisson solvers
1308 call poisson_end(psolver)
1309 call poisson_init(psolver, namespace, space, gr%der, mc, gr%stencil, qtot, verbose=.false.)
1310
1311 call kpoints_lattice_vectors_update(kpoints, new_latt)
1312
1314
1317
1318end module ion_dynamics_oct_m
1319
1320!! Local Variables:
1321!! mode: f90
1322!! coding: utf-8
1323!! 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:181
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_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
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:538
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
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
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)