29 use,
intrinsic :: iso_fortran_env
77 integer,
parameter :: &
100 logical :: relax_cell
101 logical :: constant_velocity
102 integer :: thermostat
104 real(real64) :: current_temperature
106 real(real64),
allocatable :: oldforce(:, :)
109 real(real64),
allocatable :: old_pos(:, :)
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(:)
118 real(real64) :: pressure
121 logical :: symmetrize = .false.
122 type(symmetrizer_t),
pointer :: symm
125 type(nose_hoover_t) :: nh(1:2)
126 type(tdf_t) :: temperature_function
129 logical :: drive_ions
130 type(ion_td_displacement_t),
allocatable :: td_displacements(:)
131 type(ions_t),
pointer :: ions_t0
133 real(real64),
public :: ionic_scale
143 real(real64),
allocatable :: pos(:, :)
144 real(real64),
allocatable :: vel(:, :)
145 real(real64),
allocatable :: old_pos(:, :)
146 type(nose_hoover_t) :: nh(1:2)
151 real(real64),
allocatable :: pos(:, :)
152 real(real64),
allocatable :: vel(:, :)
153 real(real64),
allocatable :: old_pos(:, :)
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
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
171 character(len=100) :: temp_function_name
172 logical :: have_velocities
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))
223 write(
message(1),
'(a)')
'Input: TDIonicTimeScale must be positive.'
241 call parse_variable(namespace,
'IonsConstantVelocity', .false., this%constant_velocity)
244 if (this%constant_velocity)
then
247 this%drive_ions = .
true.
271 if (
parse_block(namespace,
'IonsTimeDependentDisplacements', blk) == 0)
then
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.
281 this%td_displacements(iatom)%move = .
true.
282 if (iatom < 1 .and. iatom > ions%natoms)
then
287 call tdf_read(this%td_displacements(iatom)%fx, namespace, trim(expression), ierr)
289 write(
message(1),
'(3A)')
'Could not find "', trim(expression),
'" in the TDFunctions block:'
295 call tdf_read(this%td_displacements(iatom)%fy, namespace, trim(expression), ierr)
297 write(
message(1),
'(3A)')
'Could not find "', trim(expression),
'" in the TDFunctions block:'
302 call tdf_read(this%td_displacements(iatom)%fz, namespace, trim(expression), ierr)
304 write(
message(1),
'(3A)')
'Could not find "', trim(expression),
'" in the TDFunctions block:'
310 safe_allocate(this%ions_t0)
333 call parse_variable(namespace,
'Thermostat', thermo_none, this%thermostat)
337 if (this%thermostat /= thermo_none)
then
339 have_velocities = .
true.
341 if (this%drive_ions)
then
342 call messages_write(
'You cannot use a Thermostat and IonsConstantVelocity or IonsTimeDependentDisplacements')
359 call parse_variable(namespace,
'TemperatureFunction',
'temperature', temp_function_name)
361 call tdf_read(this%temperature_function, namespace, temp_function_name, ierr)
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."
381 this%nh(2)%mass = this%nh(1)%mass
386 safe_allocate(this%old_pos(1:ions%space%dim, 1:ions%natoms))
388 this%old_pos = ions%pos
410 have_velocities = .
true.
417 do i = 1, ions%natoms
420 sigma =
sqrt(temperature / ions%mass(i))
426 call mpi_world%bcast(ions%vel(:, i), ions%space%dim, mpi_double_precision, 0)
433 call ions%update_kinetic_energy()
434 kin1 = ions%kinetic_energy
436 xx = ions%center_of_mass_vel()
437 do i = 1, ions%natoms
438 ions%vel(:, i) = ions%vel(:, i) - xx
441 call ions%update_kinetic_energy()
442 kin2 = ions%kinetic_energy
445 do i = 1, ions%natoms
446 ions%vel(:, i) =
sqrt(kin1/kin2)*ions%vel(:, i)
450 call ions%update_kinetic_energy()
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> =', &
457 write(
message(3),
'(2x,a,f8.4,1x,a)')
'3/2 k_B T =', &
516 have_velocities = .
true.
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
524 do i = 1, ions%natoms
535 call ions%update_kinetic_energy()
545 call parse_variable(namespace,
'MoveIons', have_velocities, this%move_ions)
548 if (this%move_ions .and. ions%space%periodic_dim == 1)
then
550 'Moving ions for a 1D periodic system is not allowed, as forces are incorrect.')
553 if (this%ions_move())
then
554 safe_allocate(this%oldforce(1:ions%space%dim, 1:ions%natoms))
557 if (ions%space%is_periodic())
then
568 call parse_variable(namespace,
'CellDynamics', .false., this%relax_cell)
571 if (this%cell_relax())
then
573 message(1) =
"Cell dynamics not supported on GPUs."
577 periodic_dim = ions%space%periodic_dim
578 ncomp = periodic_dim * periodic_dim
579 safe_allocate(this%cell_force(1:ncomp))
581 safe_allocate(this%old_cell_force(1:ncomp))
582 this%old_cell_force =
m_zero
583 safe_allocate(this%cell_vel(1:ncomp))
586 safe_allocate(this%strain(1:ncomp))
589 do i = 1, periodic_dim
590 do j = i, periodic_dim
591 if(i == j) this%strain(ncomp) =
m_one
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)
613 this%relax_cell = .false.
626 safe_deallocate_a(this%oldforce)
628 if (this%thermostat /= thermo_none)
then
629 call tdf_end(this%temperature_function)
632 if (this%drive_ions .and.
allocated(this%td_displacements))
then
633 if (any(this%td_displacements(1:this%ions_t0%natoms)%move))
then
638 safe_deallocate_a(this%td_displacements)
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)
654 type(
ions_t),
intent(inout) :: ions
655 real(real64),
intent(in) :: time
656 real(real64),
intent(in) :: dt
665 if (this%drive_ions)
then
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))
682 if (this%ions_move())
then
687 if (this%cell_relax())
then
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))
698 ions%tot_force(:, iatom) = ions%latt%red_to_cart(ions%tot_force(:, iatom))
703 call ions%fold_atoms_into_cell()
713 real(real64),
intent(in) :: time
719 if (this%thermostat /= thermo_none)
then
722 if (this%current_temperature <
m_zero)
then
723 write(
message(1),
'(a, f10.3, 3a, f10.3, 3a)') &
724 "Negative temperature (", &
731 this%current_temperature =
m_zero
742 type(
ions_t),
intent(inout) :: ions
743 real(real64),
intent(in) :: time
744 real(real64),
intent(in) :: dt
747 real(real64) :: dr(3)
751 assert(this%drive_ions)
753 do iatom = 1, ions%natoms
754 if (ions%fixed(iatom)) cycle
756 if (this%constant_velocity)
then
757 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom)
759 else if (
allocated(this%td_displacements))
then
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) /)
766 ions%pos(:, iatom) = this%ions_t0%pos(:, iatom) + dr(1:ions%space%dim)
782 type(
ions_t),
intent(inout) :: ions
783 real(real64),
intent(in) :: time
784 real(real64),
intent(in) :: dt
790 assert(.not. this%drive_ions)
794 do iatom = 1, ions%natoms
795 if (ions%fixed(iatom)) cycle
797 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom) + &
798 m_half*dt**2 / ions%mass(iatom) * ions%tot_force(:, iatom)
800 this%oldforce(:, iatom) = ions%tot_force(:, iatom)
812 do iatom = 1, ions%natoms
813 if (ions%fixed(iatom)) cycle
815 ions%pos(:, iatom) = ions%pos(:, iatom) +
m_half*dt*ions%vel(:, iatom)
827 type(
ions_t),
intent(inout) :: ions
828 real(real64),
intent(in) :: time
829 real(real64),
intent(in) :: dt
832 integer :: idir, jdir, comp
833 real(real64) :: rlattice_change(ions%space%periodic_dim*ions%space%periodic_dim)
837 rlattice_change = dt * this%cell_vel +
m_half*dt**2 * this%cell_force
840 do idir = 1, ions%space%periodic_dim
842 this%cell_force(jdir + (idir-1)*ions%space%periodic_dim)), &
843 jdir = 1, ions%space%periodic_dim)
845 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
847 write(
message(1),
'(a,3a,a)')
' Cell vel [', &
849 do idir = 1, ions%space%periodic_dim
851 this%cell_vel(jdir+ (idir-1)*ions%space%periodic_dim)), &
852 jdir = 1, ions%space%periodic_dim)
854 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
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)
865 this%old_cell_force = this%cell_force
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)
871 call ions%update_lattice_vectors(ions%latt, this%symmetrize)
880 type(
ions_t),
intent(inout) :: ions
882 real(real64) :: g1, g2, ss, uk, dt, temp
888 call ions%update_kinetic_energy()
889 uk = ions%kinetic_energy
891 temp = this%current_temperature
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)
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
905 ions%vel = ss*ions%vel
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)
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
924 type(
ions_t),
intent(inout) :: ions
925 logical,
optional,
intent(out) :: atoms_moved
928 real(real64) :: scal, temp
931 if (this%drive_ions)
return
935 if (
present(atoms_moved)) atoms_moved = this%thermostat ==
thermo_nh
940 do iatom = 1, ions%natoms
941 if (ions%fixed(iatom)) cycle
943 ions%vel(:, iatom) = ions%vel(:, iatom) &
944 + this%dt/ions%mass(iatom) *
m_half * (this%oldforce(:, iatom) + &
945 ions%tot_force(:, iatom))
951 do iatom = 1, ions%natoms
952 if (ions%fixed(iatom)) cycle
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)
965 scal =
sqrt(this%current_temperature/temp)
966 ions%vel = scal*ions%vel
970 if (this%cell_relax())
then
971 this%cell_vel = this%cell_vel + this%dt *
m_half * (this%old_cell_force + this%cell_force)
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
994 do iatom = 1, ions%natoms
995 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
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)
1006 do iatom = 1, ions%natoms
1007 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
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
1029 do iatom = 1, ions%natoms
1030 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
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))
1041 do iatom = 1, ions%natoms
1042 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
1052 type(
ions_t),
intent(in) :: ions
1055 if (.not. this%ions_move())
return
1059 safe_allocate(state%pos(1:ions%space%dim, 1:ions%natoms))
1060 safe_allocate(state%vel(1:ions%space%dim, 1:ions%natoms))
1062 state%pos = ions%pos
1063 state%vel = ions%vel
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
1079 type(
ions_t),
intent(inout) :: ions
1082 assert(.not. this%cell_relax())
1083 if (.not. this%ions_move())
return
1087 ions%pos = state%pos
1088 ions%vel = state%vel
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)
1097 safe_deallocate_a(state%pos)
1098 safe_deallocate_a(state%vel)
1108 ions_move = this%move_ions
1118 drive_ions = this%drive_ions
1127 cell_dynamics = this%relax_cell
1136 is_active = this%relax_cell .or. this%move_ions
1144 type(ions_t),
intent(in) :: ions
1147 real(real64) :: kinetic_energy
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)
1154 temperature = m_two/m_three*kinetic_energy/ions%natoms
1163 if (this%move_ions)
then
1164 this%move_ions = .false.
1176 this%move_ions = .
true.
1182 type(restart_t),
intent(in) :: restart
1183 integer,
intent(out) :: ierr
1187 if (
allocated(this%oldforce))
then
1188 call drestart_write_binary(restart,
"ion_dynamics_oldforce",
size(this%oldforce), &
1189 this%oldforce, ierr)
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)
1203 type(restart_t),
intent(in) :: restart
1204 integer,
intent(out) :: ierr
1208 if (
allocated(this%oldforce))
then
1209 call drestart_read_binary(restart,
"ion_dynamics_oldforce",
size(this%oldforce), &
1210 this%oldforce, ierr)
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)
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
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)
1237 inv_latt = rlattice(1:space%periodic_dim, 1:space%periodic_dim)
1238 call lalg_inverse(space%periodic_dim, inv_latt,
'dir')
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
1244 cell_force = matmul(tmp_stress, transpose(inv_latt)) * rcell_volume
1247 do idir = 1, space%periodic_dim
1248 do jdir = 1, space%periodic_dim
1249 this%cell_force(comp) = cell_force(idir, jdir)
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)
1260 call messages_info(1+space%periodic_dim, namespace=global_namespace)
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
1274 real(real64) :: length(1:space%dim)
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))
1284 call box%regenerate(space%dim, new_latt%rlattice, length, namespace)
1286 call messages_not_implemented(
"Grid regeneration for non-parallelepiped boxes", namespace=namespace)
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
1305 call grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
1308 call poisson_end(psolver)
1309 call poisson_init(psolver, namespace, space, gr%der, mc, gr%stencil, qtot, verbose=.false.)
1311 call kpoints_lattice_vectors_update(kpoints, new_latt)
Functions to generate random numbers.
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
double exp(double __x) __attribute__((__nothrow__
pure logical function, public accel_is_enabled()
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module implements the underlying real-space grid.
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.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
logical function, public parse_is_defined(namespace, name)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
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)
subroutine, public tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
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