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.
280 this%td_displacements(iatom)%move = .
true.
283 call tdf_read(this%td_displacements(iatom)%fx, namespace, trim(expression), ierr)
285 write(
message(1),
'(3A)')
'Could not find "', trim(expression),
'" in the TDFunctions block:'
291 call tdf_read(this%td_displacements(iatom)%fy, namespace, trim(expression), ierr)
293 write(
message(1),
'(3A)')
'Could not find "', trim(expression),
'" in the TDFunctions block:'
298 call tdf_read(this%td_displacements(iatom)%fz, namespace, trim(expression), ierr)
300 write(
message(1),
'(3A)')
'Could not find "', trim(expression),
'" in the TDFunctions block:'
306 safe_allocate(this%ions_t0)
329 call parse_variable(namespace,
'Thermostat', thermo_none, this%thermostat)
333 if (this%thermostat /= thermo_none)
then
335 have_velocities = .
true.
337 if (this%drive_ions)
then
338 call messages_write(
'You cannot use a Thermostat and IonsConstantVelocity or IonsTimeDependentDisplacements')
355 call parse_variable(namespace,
'TemperatureFunction',
'temperature', temp_function_name)
357 call tdf_read(this%temperature_function, namespace, temp_function_name, ierr)
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."
377 this%nh(2)%mass = this%nh(1)%mass
382 safe_allocate(this%old_pos(1:ions%space%dim, 1:ions%natoms))
384 this%old_pos = ions%pos
406 have_velocities = .
true.
413 do i = 1, ions%natoms
416 sigma =
sqrt(temperature / ions%mass(i))
422 call mpi_world%bcast(ions%vel(:, i), ions%space%dim, mpi_double_precision, 0)
429 call ions%update_kinetic_energy()
430 kin1 = ions%kinetic_energy
432 xx = ions%center_of_mass_vel()
433 do i = 1, ions%natoms
434 ions%vel(:, i) = ions%vel(:, i) - xx
437 call ions%update_kinetic_energy()
438 kin2 = ions%kinetic_energy
440 do i = 1, ions%natoms
441 ions%vel(:, i) =
sqrt(kin1/kin2)*ions%vel(:, i)
444 call ions%update_kinetic_energy()
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> =', &
451 write(
message(3),
'(2x,a,f8.4,1x,a)')
'3/2 k_B T =', &
510 have_velocities = .
true.
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
518 do i = 1, ions%natoms
529 call ions%update_kinetic_energy()
539 call parse_variable(namespace,
'MoveIons', have_velocities, this%move_ions)
542 if (this%move_ions .and. ions%space%periodic_dim == 1)
then
544 'Moving ions for a 1D periodic system is not allowed, as forces are incorrect.')
547 if (this%ions_move())
then
548 safe_allocate(this%oldforce(1:ions%space%dim, 1:ions%natoms))
551 if (ions%space%is_periodic())
then
562 call parse_variable(namespace,
'CellDynamics', .false., this%relax_cell)
565 if (this%cell_relax())
then
567 message(1) =
"Cell dynamics not supported on GPUs."
571 periodic_dim = ions%space%periodic_dim
572 ncomp = periodic_dim * periodic_dim
573 safe_allocate(this%cell_force(1:ncomp))
575 safe_allocate(this%old_cell_force(1:ncomp))
576 this%old_cell_force =
m_zero
577 safe_allocate(this%cell_vel(1:ncomp))
580 safe_allocate(this%strain(1:ncomp))
583 do i = 1, periodic_dim
584 do j = i, periodic_dim
585 if(i == j) this%strain(ncomp) =
m_one
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)
607 this%relax_cell = .false.
620 safe_deallocate_a(this%oldforce)
622 if (this%thermostat /= thermo_none)
then
623 call tdf_end(this%temperature_function)
626 if (this%drive_ions .and.
allocated(this%td_displacements))
then
627 if (any(this%td_displacements(1:this%ions_t0%natoms)%move))
then
632 safe_deallocate_a(this%td_displacements)
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)
648 type(
ions_t),
intent(inout) :: ions
649 real(real64),
intent(in) :: time
650 real(real64),
intent(in) :: dt
659 if (this%drive_ions)
then
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))
676 if (this%ions_move())
then
681 if (this%cell_relax())
then
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))
692 ions%tot_force(:, iatom) = ions%latt%red_to_cart(ions%tot_force(:, iatom))
697 call ions%fold_atoms_into_cell()
707 real(real64),
intent(in) :: time
713 if (this%thermostat /= thermo_none)
then
716 if (this%current_temperature <
m_zero)
then
717 write(
message(1),
'(a, f10.3, 3a, f10.3, 3a)') &
718 "Negative temperature (", &
725 this%current_temperature =
m_zero
736 type(
ions_t),
intent(inout) :: ions
737 real(real64),
intent(in) :: time
738 real(real64),
intent(in) :: dt
741 real(real64) :: dr(3)
745 assert(this%drive_ions)
747 do iatom = 1, ions%natoms
748 if (ions%fixed(iatom)) cycle
750 if (this%constant_velocity)
then
751 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom)
753 else if (
allocated(this%td_displacements))
then
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) /)
760 ions%pos(:, iatom) = this%ions_t0%pos(:, iatom) + dr(1:ions%space%dim)
776 type(
ions_t),
intent(inout) :: ions
777 real(real64),
intent(in) :: time
778 real(real64),
intent(in) :: dt
784 assert(.not. this%drive_ions)
788 do iatom = 1, ions%natoms
789 if (ions%fixed(iatom)) cycle
791 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom) + &
792 m_half*dt**2 / ions%mass(iatom) * ions%tot_force(:, iatom)
794 this%oldforce(:, iatom) = ions%tot_force(:, iatom)
806 do iatom = 1, ions%natoms
807 if (ions%fixed(iatom)) cycle
809 ions%pos(:, iatom) = ions%pos(:, iatom) +
m_half*dt*ions%vel(:, iatom)
821 type(
ions_t),
intent(inout) :: ions
822 real(real64),
intent(in) :: time
823 real(real64),
intent(in) :: dt
826 integer :: idir, jdir, comp
827 real(real64) :: rlattice_change(ions%space%periodic_dim*ions%space%periodic_dim)
831 rlattice_change = dt * this%cell_vel +
m_half*dt**2 * this%cell_force
834 do idir = 1, ions%space%periodic_dim
836 this%cell_force(jdir + (idir-1)*ions%space%periodic_dim)), &
837 jdir = 1, ions%space%periodic_dim)
839 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
841 write(
message(1),
'(a,3a,a)')
' Cell vel [', &
843 do idir = 1, ions%space%periodic_dim
845 this%cell_vel(jdir+ (idir-1)*ions%space%periodic_dim)), &
846 jdir = 1, ions%space%periodic_dim)
848 call messages_info(1+ions%space%periodic_dim, namespace=ions%namespace)
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)
859 this%old_cell_force = this%cell_force
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)
865 call ions%update_lattice_vectors(ions%latt, this%symmetrize)
874 type(
ions_t),
intent(inout) :: ions
876 real(real64) :: g1, g2, ss, uk, dt, temp
882 call ions%update_kinetic_energy()
883 uk = ions%kinetic_energy
885 temp = this%current_temperature
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)
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
899 ions%vel = ss*ions%vel
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)
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
918 type(
ions_t),
intent(inout) :: ions
919 logical,
optional,
intent(out) :: atoms_moved
925 if (this%drive_ions)
return
929 if (
present(atoms_moved)) atoms_moved = this%thermostat ==
thermo_nh
934 do iatom = 1, ions%natoms
935 if (ions%fixed(iatom)) cycle
937 ions%vel(:, iatom) = ions%vel(:, iatom) &
938 + this%dt/ions%mass(iatom) *
m_half * (this%oldforce(:, iatom) + &
939 ions%tot_force(:, iatom))
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)
957 ions%vel = scal*ions%vel
960 if (this%cell_relax())
then
961 this%cell_vel = this%cell_vel + this%dt *
m_half * (this%old_cell_force + this%cell_force)
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
982 do iatom = 1, ions%natoms
983 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
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)
994 do iatom = 1, ions%natoms
995 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
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
1017 do iatom = 1, ions%natoms
1018 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
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))
1029 do iatom = 1, ions%natoms
1030 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
1040 type(
ions_t),
intent(in) :: ions
1043 if (.not. this%ions_move())
return
1047 safe_allocate(state%pos(1:ions%space%dim, 1:ions%natoms))
1048 safe_allocate(state%vel(1:ions%space%dim, 1:ions%natoms))
1050 state%pos = ions%pos
1051 state%vel = ions%vel
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
1067 type(
ions_t),
intent(inout) :: ions
1070 assert(.not. this%cell_relax())
1071 if (.not. this%ions_move())
return
1075 ions%pos = state%pos
1076 ions%vel = state%vel
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)
1085 safe_deallocate_a(state%pos)
1086 safe_deallocate_a(state%vel)
1096 ions_move = this%move_ions
1106 drive_ions = this%drive_ions
1115 cell_dynamics = this%relax_cell
1124 is_active = this%relax_cell .or. this%move_ions
1132 type(ions_t),
intent(in) :: ions
1135 real(real64) :: kinetic_energy
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)
1142 temperature = m_two/m_three*kinetic_energy/ions%natoms
1151 if (this%move_ions)
then
1152 this%move_ions = .false.
1164 this%move_ions = .
true.
1170 type(restart_t),
intent(in) :: restart
1171 integer,
intent(out) :: ierr
1175 if (
allocated(this%oldforce))
then
1176 call drestart_write_binary(restart,
"ion_dynamics_oldforce",
size(this%oldforce), &
1177 this%oldforce, ierr)
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)
1191 type(restart_t),
intent(in) :: restart
1192 integer,
intent(out) :: ierr
1196 if (
allocated(this%oldforce))
then
1197 call drestart_read_binary(restart,
"ion_dynamics_oldforce",
size(this%oldforce), &
1198 this%oldforce, ierr)
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)
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
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)
1225 inv_latt = rlattice(1:space%periodic_dim, 1:space%periodic_dim)
1226 call lalg_inverse(space%periodic_dim, inv_latt,
'dir')
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
1232 cell_force = matmul(tmp_stress, transpose(inv_latt)) * rcell_volume
1235 do idir = 1, space%periodic_dim
1236 do jdir = 1, space%periodic_dim
1237 this%cell_force(comp) = cell_force(idir, jdir)
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)
1248 call messages_info(1+space%periodic_dim, namespace=global_namespace)
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
1262 real(real64) :: length(1:space%dim)
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))
1272 call box%regenerate(space%dim, new_latt%rlattice, length, namespace)
1274 call messages_not_implemented(
"Grid regeneration for non-parallelepiped boxes", namespace=namespace)
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
1293 call grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
1296 call poisson_end(psolver)
1297 call poisson_init(psolver, namespace, space, gr%der, mc, gr%stencil, qtot, verbose=.false.)
1299 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_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