37 use,
intrinsic :: iso_fortran_env
71 type(lattice_vectors_t) :: latt
74 type(atom_t),
allocatable :: atom(:)
76 type(symmetries_t) :: symm
78 type(distributed_t) :: atoms_dist
80 integer,
allocatable :: map_symm_atoms(:,:)
81 integer,
allocatable :: inv_map_symm_atoms(:,:)
83 real(real64),
allocatable :: equilibrium_pos(:,:)
89 real(real64),
allocatable :: pos_displacements(:,:)
90 real(real64),
allocatable :: vel_displacements(:,:)
94 class(species_wrapper_t),
allocatable :: species(:)
95 logical :: only_user_def
96 logical,
private :: species_time_dependent
98 logical :: force_total_enforce
99 type(ion_interaction_t) :: ion_interaction
102 logical,
private :: apply_global_force
103 type(tdf_t),
private :: global_force_function
106 generic ::
assignment(=) => copy
139 procedure ions_constructor
146 type(namespace_t),
intent(in) :: namespace
147 logical,
optional,
intent(in) :: print_info
148 type(lattice_vectors_t),
optional,
intent(out) :: latt_inp
149 class(ions_t),
pointer :: ions
151 type(read_coords_info) :: xyz
152 integer :: ia, ierr, idir
153 character(len=100) :: function_name
154 real(real64) :: mindist
155 real(real64),
allocatable :: factor(:)
156 integer,
allocatable :: site_type(:)
157 logical,
allocatable :: spherical_site(:)
158 real(real64),
parameter :: threshold = 1e-5_real64
159 type(species_factory_t) :: factory
161 integer :: displacement_mode
162 real(real64) :: amp_pos, amp_vel
164 logical :: in_ensemble
170 ions%namespace = namespace
172 ions%space =
space_t(namespace)
185 message(1) =
"Coordinates have not been defined."
194 safe_allocate(ions%atom(1:ions%natoms))
195 do ia = 1, ions%natoms
196 call atom_init(ions%atom(ia), ions%space%dim, xyz%atom(ia)%label)
197 ions%pos(:,ia) = xyz%atom(ia)%x(1:ions%space%dim)
199 ions%fixed(ia) = .not. xyz%atom(ia)%move
203 if (
allocated(xyz%latvec))
then
213 do ia = 1, ions%natoms
214 ions%pos(:, ia) = ions%latt%red_to_cart(ions%pos(:, ia))
221 safe_allocate_source_a(ions%equilibrium_pos, ions%pos)
226 do ia = 1, ions%natoms
227 ions%mass(ia) = ions%atom(ia)%species%get_mass()
228 ions%charge(ia) = ions%atom(ia)%species%get_zval()
242 call parse_variable(namespace,
'InitialDisplacementMode', 0, displacement_mode)
252 call parse_variable(namespace,
'InitialDisplacementAmplitudePos', 1.0_real64, amp_pos)
262 call parse_variable(namespace,
'InitialDisplacementAmplitudeVel', 1.0_real64, amp_vel)
264 if (displacement_mode>0)
then
266 message(1) =
"Random displacements are disabled when InitialDisplacementMode > 0."
269 safe_allocate(ions%pos_displacements(1:ions%space%dim, 1:ions%natoms))
270 safe_allocate(ions%vel_displacements(1:ions%space%dim, 1:ions%natoms))
272 call ions%single_mode_displacements(displacement_mode, amp_pos, amp_vel)
275 ions%pos = ions%pos + ions%pos_displacements
283 if(in_ensemble .and. displacement_mode==0)
then
292 call parse_variable(namespace,
'EnsembleTemperature', 0.0_real64, t)
294 safe_allocate(ions%pos_displacements(1:ions%space%dim, 1:ions%natoms))
295 safe_allocate(ions%vel_displacements(1:ions%space%dim, 1:ions%natoms))
296 call ions%init_random_displacements(t)
299 ions%pos = ions%pos + ions%pos_displacements
310 if (
present(latt_inp))
then
317 if (ions%space%has_mixed_periodicity())
then
318 safe_allocate(factor(ions%space%dim))
319 do idir = 1, ions%space%periodic_dim
322 do idir = ions%space%periodic_dim + 1, ions%space%dim
323 factor(idir) =
m_one/norm2(ions%latt%rlattice(1:ions%space%dim, idir))
325 call ions%latt%scale(factor)
326 safe_deallocate_a(factor)
330 if (ions%natoms > 1)
then
332 if (mindist < threshold)
then
333 write(
message(1),
'(a)')
"Some of the atoms seem to sit too close to each other."
334 write(
message(2),
'(a)')
"Please review your input files and the output geometry (in 'static/')."
335 write(
message(3),
'(a, f12.6, 1x, a)')
"Minimum distance = ", &
341 call ions%write_xyz(trim(
static_dir)//
'/geometry')
345 message(1) =
"It cannot be correct to run with physical atoms so close."
351 safe_allocate(spherical_site(1:ions%natoms))
352 safe_allocate(site_type(1:ions%natoms))
353 do ia = 1, ions%natoms
354 select type(spec => ions%atom(ia)%species)
356 spherical_site(ia) = .false.
358 spherical_site(ia) = .false.
360 spherical_site(ia) = .false.
362 spherical_site(ia) = .false.
364 spherical_site(ia) = .
true.
367 site_type(ia) = ions%atom(ia)%species%get_index()
370 ions%symm =
symmetries_t(ions%namespace, ions%space, ions%latt, ions%natoms, ions%pos, site_type, spherical_site)
372 safe_deallocate_a(spherical_site)
373 safe_deallocate_a(site_type)
376 call ions%generate_mapping_symmetric_atoms()
388 call parse_variable(namespace,
'ForceTotalEnforce', .false., ions%force_total_enforce)
404 ions%apply_global_force = .
true.
406 call parse_variable(namespace,
'TDGlobalForce',
'none', function_name)
407 call tdf_read(ions%global_force_function, namespace, trim(function_name), ierr)
410 call messages_write(
"You have enabled the GlobalForce option but Octopus could not find")
411 call messages_write(
"the '"//trim(function_name)//
"' function in the TDFunctions block.")
417 ions%apply_global_force = .false.
428 type(
ions_t),
intent(inout) :: ions
430 logical,
optional,
intent(in) :: print_info
432 logical :: print_info_, spec_user_defined
433 integer :: i, j, k, ispin
439 if (
present(print_info))
then
440 print_info_ = print_info
444 atoms1:
do i = 1, ions%natoms
448 ions%nspecies = ions%nspecies + 1
452 allocate(ions%species(1:ions%nspecies))
456 ions%only_user_def = .
true.
457 atoms2:
do i = 1, ions%natoms
462 ions%species(k)%s => factory%create_from_input(ions%namespace, ions%atom(j)%get_label(), k)
464 select type(spec => ions%species(k)%s)
468 ions%only_user_def = .false.
471 select type(spec => ions%species(k)%s)
473 if (ions%space%dim /= 3)
then
474 message(1) =
"Pseudopotentials may only be used with Dimensions = 3."
479 if (ions%space%is_periodic() .and. ions%space%periodic_dim /= 2)
then
480 message(1) =
"Periodic jelium slab can only be used if PeriodicDim = 2"
491 ispin = min(2, ispin)
493 if (print_info_)
then
496 do i = 1, ions%nspecies
497 spec => ions%species(i)%s
498 call spec%build(ions%namespace, ispin, ions%space%dim, print_info=print_info_)
500 if (print_info_)
then
513 call parse_variable(ions%namespace,
'SpeciesTimeDependent', .false., ions%species_time_dependent)
515 do i = 1,ions%nspecies
516 select type(spec=>ions%species(i)%s)
518 spec_user_defined = .
true.
521 if (ions%species_time_dependent .and. .not. spec_user_defined)
then
526 do i = 1, ions%natoms
527 do j = 1, ions%nspecies
540 class(
ions_t),
intent(out) :: ions_out
541 class(
ions_t),
intent(in) :: ions_in
547 ions_out%latt = ions_in%latt
549 ions_out%natoms = ions_in%natoms
550 safe_allocate(ions_out%atom(1:ions_out%natoms))
551 ions_out%atom = ions_in%atom
553 ions_out%nspecies = ions_in%nspecies
554 allocate(ions_out%species(1:ions_out%nspecies))
555 ions_out%species = ions_in%species
557 ions_out%only_user_def = ions_in%only_user_def
561 safe_allocate(ions_out%map_symm_atoms(1:ions_in%natoms, 1:ions_in%symm%nops))
562 ions_out%map_symm_atoms = ions_in%map_symm_atoms
563 safe_allocate(ions_out%inv_map_symm_atoms(1:ions_in%natoms, 1:ions_in%symm%nops))
564 ions_out%inv_map_symm_atoms = ions_in%inv_map_symm_atoms
572 class(
ions_t),
intent(inout) :: this
586 class(
ions_t),
target,
intent(inout) :: this
591 select type (interaction)
609 class(
ions_t),
intent(inout) :: ions
610 real(real64),
intent(in) :: T
613 integer(int64) :: seed
618 real(real64),
allocatable :: sigmaP(:)
619 real(real64),
allocatable :: muP(:)
620 real(real64),
allocatable :: sigmaQ(:)
621 real(real64),
allocatable :: muQ(:)
622 real(real64),
allocatable :: genP(:)
623 real(real64),
allocatable :: genQ(:)
626 real(real64) :: beta_half
627 real(real64),
allocatable :: pos_displacements(:), vel_displacements(:)
628 real(real64),
allocatable :: l_m(:), prefactor(:)
630 integer :: imode, idim, iatom, i, iline
631 integer :: num_real_modes
633 real(real64),
parameter :: low_temperature_tolerance = 1.0e-6_real64
639 seed = ions%namespace%get_hash32()
641 message(1) =
"Create initial random displacements for ions."
644 write(
message(1),
'("namespace = ",A)') ions%namespace%get()
645 write(
message(2),
'("seed = ",I0)') seed
650 call phonons%init(ions%namespace, ions%space%dim, ions%natoms, ions%space%periodic_dim > 0)
652 num_real_modes = phonons%num_modes
654 if (num_real_modes == 0)
then
656 ions%pos_displacements =
m_zero
657 ions%vel_displacements =
m_zero
665 call wigner%init(num_real_modes, seed)
668 safe_allocate(sigmap(1:num_real_modes))
669 safe_allocate(sigmaq(1:num_real_modes))
670 safe_allocate(mup(1:num_real_modes))
671 safe_allocate(muq(1:num_real_modes))
672 safe_allocate(genp(1:num_real_modes))
673 safe_allocate(genq(1:num_real_modes))
675 safe_allocate(pos_displacements(1:(ions%space%dim*ions%natoms)))
676 safe_allocate(vel_displacements(1:(ions%space%dim*ions%natoms)))
677 safe_allocate(prefactor(1:(ions%space%dim*ions%natoms)))
679 safe_allocate(l_m(1:num_real_modes))
682 if (t < low_temperature_tolerance)
then
684 sigmap = 1.0_real64/2.0_real64
685 sigmaq = 1.0_real64/2.0_real64
688 sigmap = 1.0_real64/
sqrt((2.0_real64 *
tanh(beta_half * phonons%frequencies)))
689 sigmaq = 1.0_real64/
sqrt((2.0_real64 *
tanh(beta_half * phonons%frequencies)))
694 l_m =
sqrt(2.0_real64/(
unit_amu%factor * phonons%frequencies)) * phonons%alpha
697 do iatom = 1, ions%natoms
698 do idim = 1, ions%space%dim
699 prefactor(i) = 1.0_real64/
sqrt(ions%mass(iatom)/
unit_amu%factor * phonons%num_super)
705 genq = wigner%get(sigmaq, muq,
wigner_q) * l_m
710 ions%pos_displacements =
m_zero
711 ions%vel_displacements =
m_zero
715 do imode = 1, num_real_modes
717 pos_displacements = prefactor(:) * phonons%eigenvectors(:, imode) * genq(imode)
718 vel_displacements = prefactor(:) * phonons%eigenvectors(:, imode) * genp(imode)
721 write(
message(iline),
'("Displacements for mode ",I4)') imode
722 do iatom=1, ions%natoms
724 write(
message(iline),
'(3E15.5)') pos_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
727 write(
message(iline),
'("Velocities for mode ",I4)') imode
728 do iatom=1, ions%natoms
730 write(
message(iline),
'(3E15.5)') vel_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
735 call blas_axpy(phonons%dim, 1.0_real64, pos_displacements(1), 1, ions%pos_displacements(1,1), 1)
736 call blas_axpy(phonons%dim, 1.0_real64, vel_displacements(1), 1, ions%vel_displacements(1,1), 1)
740 safe_deallocate_a(sigmap)
741 safe_deallocate_a(sigmaq)
742 safe_deallocate_a(mup)
743 safe_deallocate_a(muq)
744 safe_deallocate_a(genp)
745 safe_deallocate_a(genq)
746 safe_deallocate_a(l_m)
748 safe_deallocate_a(prefactor)
749 safe_deallocate_a(pos_displacements)
750 safe_deallocate_a(vel_displacements)
762 class(
ions_t),
intent(inout) :: ions
763 integer,
intent(in) :: mode
764 real(real64),
optional,
intent(in) :: amplitude_pos
765 real(real64),
optional,
intent(in) :: amplitude_vel
768 integer :: num_real_modes, i, iatom, idim
769 real(real64) :: l_m, genP, genQ
771 real(real64),
allocatable :: pos_displacements(:), vel_displacements(:), prefactor(:)
776 write(
message(1),
'("Create displacements for mode ", I4)') mode
780 call phonons%init(ions%namespace, ions%space%dim, ions%natoms, ions%space%periodic_dim > 0)
782 num_real_modes = phonons%num_modes
784 if (mode > num_real_modes)
then
785 write(
message(1),
'("Requested mode ",I0," exceeds number of available modes (",I0,")")') mode, num_real_modes
792 ions%pos_displacements =
m_zero
793 ions%vel_displacements =
m_zero
795 l_m =
sqrt(2.0_real64/(
unit_amu%factor * phonons%frequencies(mode))) * phonons%alpha(mode)
797 safe_allocate(pos_displacements(1:(ions%space%dim*ions%natoms)))
798 safe_allocate(vel_displacements(1:(ions%space%dim*ions%natoms)))
799 safe_allocate(prefactor(1:(ions%space%dim*ions%natoms)))
802 do iatom = 1, ions%natoms
803 do idim=1, ions%space%dim
804 prefactor(i) = 1.0_real64/
sqrt(ions%mass(iatom)/
unit_amu%factor * phonons%num_super)
810 genq = amplitude_pos * l_m
811 genp = amplitude_vel / (l_m *
unit_amu%factor)
813 pos_displacements = prefactor(:) * phonons%eigenvectors(:, mode) * genq
814 vel_displacements = prefactor(:) * phonons%eigenvectors(:, mode) * genp
816 write(
message(1),
'("Displacements for mode ",I4)') mode
818 do iatom=1, ions%natoms
819 write(
message(1),
'(3E15.5)') pos_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
822 write(
message(1),
'("Velocities for mode ",I4)') mode
824 do iatom=1, ions%natoms
825 write(
message(1),
'(3E15.5)') vel_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
830 call blas_axpy(phonons%dim, 1.0_real64, pos_displacements(1), 1, ions%pos_displacements(1,1), 1)
831 call blas_axpy(phonons%dim, 1.0_real64, vel_displacements(1), 1, ions%vel_displacements(1,1), 1)
833 safe_deallocate_a(pos_displacements)
834 safe_deallocate_a(vel_displacements)
835 safe_deallocate_a(prefactor)
843 class(
ions_t),
intent(inout) :: this
848 if(
allocated(this%vel_displacements))
then
849 this%vel = this%vel + this%vel_displacements
857 class(
ions_t),
intent(inout) :: this
858 character(len=*),
intent(in) :: label
873 class(
ions_t),
intent(in) :: partner
878 select type (interaction)
888 class(
ions_t),
intent(inout) :: partner
893 select type (interaction)
904 class(
ions_t),
intent(inout) :: this
910 do iatom = 1, this%natoms
911 this%pos(:, iatom) = this%latt%fold_into_cell(this%pos(:, iatom))
918 real(real64) function ions_min_distance(this, real_atoms_only) result(rmin)
919 class(
ions_t),
intent(in) :: this
920 logical,
optional,
intent(in) :: real_atoms_only
922 integer :: iatom, jatom, idir
923 real(real64) :: xx(this%space%dim)
924 logical :: real_atoms_only_
927 if (this%natoms == 1 .and. .not. this%space%is_periodic())
then
932 push_sub(ions_min_distance)
941 do iatom = 1, this%natoms
945 if (real_atoms_only_) cycle
947 do jatom = iatom + 1, this%natoms
951 if (real_atoms_only_) cycle
953 xx = abs(this%pos(:, iatom) - this%pos(:, jatom))
954 if (this%space%is_periodic())
then
955 xx = this%latt%cart_to_red(xx)
956 do idir = 1, this%space%periodic_dim
959 xx = this%latt%red_to_cart(xx)
961 rmin = min(norm2(xx), rmin)
965 if (.not. (this%only_user_def .and. real_atoms_only_))
then
967 do idir = 1, this%space%periodic_dim
968 rmin = min(rmin, norm2(this%latt%rlattice(:,idir)))
973 if(rmin < huge(rmin)/1e6_real64)
then
974 rmin = anint(rmin*1e6_real64)*1.0e-6_real64
977 pop_sub(ions_min_distance)
986 time_dependent = this%species_time_dependent
992 real(real64) function ions_val_charge(this, mask) result(val_charge)
993 class(
ions_t),
intent(in) :: this
994 logical,
optional,
intent(in) :: mask(:)
998 if (
present(mask))
then
999 val_charge = -sum(this%charge, mask=mask)
1001 val_charge = -sum(this%charge)
1009 class(
ions_t),
intent(in) :: this
1010 logical,
optional,
intent(in) :: mask(:)
1011 real(real64) :: dipole(this%space%dim)
1018 do ia = 1, this%natoms
1019 if (
present(mask))
then
1020 if (.not. mask(ia)) cycle
1022 dipole = dipole + this%charge(ia)*this%pos(:, ia)
1024 dipole = p_proton_charge*dipole
1031 class(
ions_t),
intent(inout) :: this
1032 real(real64),
intent(in) :: xx(this%space%dim)
1038 do iatom = 1, this%natoms
1039 this%pos(:, iatom) = this%pos(:, iatom) - xx
1047 class(
ions_t),
intent(inout) :: this
1048 real(real64),
intent(in) :: from(this%space%dim)
1049 real(real64),
intent(in) :: from2(this%space%dim)
1050 real(real64),
intent(in) :: to(this%space%dim)
1053 real(real64) :: m1(3, 3), m2(3, 3)
1054 real(real64) :: m3(3, 3), f2(3), per(3)
1055 real(real64) :: alpha, r
1059 if (this%space%dim /= 3)
then
1060 call messages_not_implemented(
"ions_rotate in other than 3 dimensions", namespace=this%namespace)
1064 m1 = diagonal_matrix(3, m_one)
1067 if (abs(to(2)) > 1d-150)
then
1068 alpha =
atan2(to(2), to(1))
1069 call rotate(m1, alpha, 3)
1071 alpha =
atan2(norm2(to(1:2)), to(3))
1072 call rotate(m1, -alpha, 2)
1075 f2 = matmul(m1, from)
1080 if (r > m_zero)
then
1087 m2 = diagonal_matrix(3, m_one)
1088 alpha =
atan2(per(1), per(2))
1089 call rotate(m2, -alpha, 3)
1092 m3 = diagonal_matrix(3, m_one)
1093 alpha =
acos(min(sum(from*to), m_one))
1094 call rotate(m3, -alpha, 2)
1097 m2 = matmul(transpose(m2), matmul(m3, m2))
1100 per = matmul(m2, matmul(m1, from2))
1101 alpha =
atan2(per(1), per(2))
1102 call rotate(m2, -alpha, 3)
1105 m1 = matmul(transpose(m1), matmul(m2, m1))
1109 do iatom = 1, this%natoms
1110 f2 = this%pos(:, iatom)
1111 this%pos(:, iatom) = matmul(m1, f2)
1118 subroutine rotate(m, angle, dir)
1119 real(real64),
intent(inout) :: m(3, 3)
1120 real(real64),
intent(in) :: angle
1121 integer,
intent(in) :: dir
1123 real(real64) :: aux(3, 3), ca, sa
1161 class(
ions_t),
intent(in) :: this
1162 real(real64),
intent(in) :: time
1163 real(real64) :: force(this%space%dim)
1169 if (this%apply_global_force)
then
1170 force(1) = tdf(this%global_force_function, time)
1177 subroutine ions_write_xyz(this, fname, append, comment, reduce_coordinates)
1178 class(
ions_t),
intent(in) :: this
1179 character(len=*),
intent(in) :: fname
1180 logical,
optional,
intent(in) :: append
1181 character(len=*),
optional,
intent(in) :: comment
1182 logical,
optional,
intent(in) :: reduce_coordinates
1184 integer :: iatom, idim, iunit
1185 character(len=6) position
1186 character(len=19) :: frmt
1187 real(real64) :: red(this%space%dim)
1189 if (.not. mpi_world%is_root())
return
1194 if (
present(append))
then
1195 if (append) position =
'append'
1197 if(.not.optional_default(reduce_coordinates, .false.))
then
1198 iunit = io_open(trim(fname)//
'.xyz', this%namespace, action=
'write', position=position)
1200 iunit = io_open(trim(fname)//
'.xyz_red', this%namespace, action=
'write', position=position)
1203 write(iunit,
'(i4)') this%natoms
1204 if (
present(comment))
then
1205 write(iunit,
'(1x,a)') comment
1207 write(iunit,
'(1x,a,a)')
'units: ', trim(units_abbrev(units_out%length_xyz_file))
1210 write(unit=frmt, fmt=
"(a5,i2.2,a4,i2.2,a6)")
"(6x,a", label_len,
",2x,", this%space%dim,
"f15.9)"
1211 do iatom = 1, this%natoms
1212 if(.not.optional_default(reduce_coordinates, .false.))
then
1213 write(unit=iunit, fmt=frmt) this%atom(iatom)%label, &
1214 (units_from_atomic(units_out%length_xyz_file, this%pos(idim, iatom)), idim=1, this%space%dim)
1216 red = this%latt%cart_to_red(this%pos(:, iatom))
1217 write(unit=iunit, fmt=frmt) this%atom(iatom)%label, (red(idim), idim=1, this%space%dim)
1220 call io_close(iunit)
1231 class(
ions_t),
intent(in) :: this
1232 character(len=*),
optional,
intent(in) :: fname
1233 character(len=*),
optional,
intent(in) :: comment
1235 integer :: iatom, idim, iunit
1236 character(len=:),
allocatable :: fname_, comment_
1237 character(len=10) :: format_string
1241 if (.not. mpi_world%is_root())
then
1246 comment_ = optional_default(comment,
"")
1247 fname_ = optional_default(fname,
"POSCAR")
1249 iunit = io_open(trim(fname_), this%namespace, action=
'write')
1251 write(iunit,
'(A)') comment_
1252 write(iunit,
'("1.0")')
1254 write(format_string,
'("(",I1,A,")")') this%space%dim,
'F15.10'
1255 do idim=1, this%space%dim
1256 write(iunit, format_string) units_from_atomic(unit_angstrom, this%latt%rlattice(:,idim))
1259 do iatom = 1, this%natoms
1260 write(iunit,
'(A)', advance=
'NO') trim(this%atom(iatom)%label)//
" "
1262 write(iunit,
'("")')
1263 write(iunit,
'(A)') repeat(
"1 ", this%natoms)
1264 write(iunit,
'("direct")')
1265 write(format_string,
'("(",I1,A,")")') this%space%dim,
'F15.10'
1266 do iatom=1, this%natoms
1267 write(iunit, format_string) this%latt%cart_to_red(this%pos(:, iatom))
1270 call io_close(iunit)
1277 use iso_fortran_env,
only: real64
1278 class(
ions_t),
intent(inout) :: this
1279 character(len=*),
intent(in) :: fname
1280 character(len=*),
optional,
intent(in) :: comment
1282 integer :: iatom, idir, iunit, ios
1283 character(len=256) :: line
1284 character(len=LABEL_LEN) :: label
1285 character(len=:),
allocatable :: coordstr
1286 real(real64) :: tmp(this%space%dim)
1290 iunit = io_open(trim(fname)//
'.xyz', this%namespace, action=
'read', position=
'rewind')
1293 read(iunit,
'(i4)') this%natoms
1294 if (
present(comment))
then
1301 do iatom = 1, this%natoms
1304 read(iunit,
'(A)', iostat=ios) line
1306 call io_close(iunit)
1307 message(1) =
"Error reading XYZ atom line"
1308 call messages_fatal(1, namespace=this%namespace)
1312 if (len_trim(line) < label_len)
then
1313 call io_close(iunit)
1314 message(1) =
"XYZ file: atom line too short for label"
1315 call messages_fatal(1, namespace=this%namespace)
1317 label = line(1:label_len)
1320 if (len(line) > label_len)
then
1321 coordstr = adjustl(line(label_len+1:))
1323 call io_close(iunit)
1324 message(1) =
"XYZ file: no coordinate data after label"
1325 call messages_fatal(1, namespace=this%namespace)
1329 read(coordstr, *, iostat=ios) (tmp(idir), idir=1, this%space%dim)
1331 call io_close(iunit)
1332 message(1) =
"XYZ file: malformed coordinate fields"
1333 call messages_fatal(1, namespace=this%namespace)
1337 this%pos(:, iatom) = units_to_atomic(units_out%length_xyz_file, tmp)
1341 call io_close(iunit)
1350 class(
ions_t),
intent(in) :: this
1351 character(len=*),
intent(in) :: dir
1353 type(lattice_iterator_t) :: latt_iter
1354 real(real64) :: radius, pos(this%space%dim)
1355 integer :: iatom, icopy, iunit
1359 radius = maxval(m_half*norm2(this%latt%rlattice(:,1:this%space%periodic_dim), dim=1))*(m_one + m_epsilon)
1360 latt_iter = lattice_iterator_t(this%latt, radius)
1362 if (mpi_world%is_root())
then
1364 iunit = io_open(trim(dir)//
'/crystal.xyz', this%namespace, action=
'write')
1366 write(iunit,
'(i9)') this%natoms*latt_iter%n_cells
1367 write(iunit,
'(a)')
'#generated by Octopus'
1369 do iatom = 1, this%natoms
1370 do icopy = 1, latt_iter%n_cells
1371 pos = units_from_atomic(units_out%length, this%pos(:, iatom) + latt_iter%get(icopy))
1372 write(iunit,
'(a, 99f12.6)') this%atom(iatom)%label, pos
1376 call io_close(iunit)
1384 class(
ions_t),
intent(in) :: this
1385 character(len=*),
intent(in) :: dir, fname
1387 integer :: iunit, iatom, idir
1388 real(real64) :: force(this%space%dim), center(this%space%dim)
1389 character(len=20) frmt
1391 if (.not. mpi_world%is_root())
return
1395 call io_mkdir(dir, this%namespace)
1396 iunit = io_open(trim(dir)//
'/'//trim(fname)//
'.bild', this%namespace, action=
'write', &
1399 write(frmt,
'(a,i0,a)')
'(a,2(', this%space%dim,
'f16.6,1x))'
1401 write(iunit,
'(a)')
'.comment : force vectors in ['//trim(units_abbrev(units_out%force))//
']'
1403 write(iunit,
'(a)')
'.color red'
1405 do iatom = 1, this%natoms
1406 center = units_from_atomic(units_out%length, this%pos(:, iatom))
1407 force = units_from_atomic(units_out%force, this%tot_force(:, iatom))
1408 write(iunit,
'(a,1x,i4,1x,a2,1x,a6,1x,f10.6,a)')
'.comment :', iatom, trim(this%atom(iatom)%label), &
1409 'force:', norm2(force),
'['//trim(units_abbrev(units_out%force))//
']'
1410 write(iunit,fmt=trim(frmt))
'.arrow', (center(idir), idir = 1, this%space%dim), &
1411 (center(idir) + force(idir), idir = 1, this%space%dim)
1415 call io_close(iunit)
1422 class(
ions_t),
intent(in) :: this
1423 character(len=*),
intent(in) :: filename
1424 logical,
optional,
intent(in) :: ascii
1426 integer :: iunit, iatom, ierr
1428 real(real64),
allocatable :: data(:, :)
1429 integer,
allocatable :: idata(:, :)
1430 character(len=MAX_PATH_LEN) :: fullname
1434 assert(this%space%dim == 3)
1436 ascii_ = optional_default(ascii, .
true.)
1438 fullname = trim(filename)//
".vtk"
1440 iunit = io_open(trim(fullname), this%namespace, action=
'write')
1442 write(iunit,
'(1a)')
'# vtk DataFile Version 3.0 '
1443 write(iunit,
'(6a)')
'Generated by octopus ', trim(conf%version),
' - git: ', &
1444 trim(conf%git_commit),
" configuration: ", trim(conf%config_time)
1447 write(iunit,
'(1a)')
'ASCII'
1449 write(iunit,
'(1a)')
'BINARY'
1452 write(iunit,
'(1a)')
'DATASET POLYDATA'
1454 write(iunit,
'(a,i9,a)')
'POINTS ', this%natoms,
' double'
1457 do iatom = 1, this%natoms
1458 write(iunit,
'(3f12.6)') this%pos(1:3, iatom)
1461 call io_close(iunit)
1462 safe_allocate(
data(1:3, 1:this%natoms))
1463 do iatom = 1, this%natoms
1464 data(1:3, iatom) = this%pos(1:3, iatom)
1466 call io_binary_write(io_workpath(fullname, this%namespace), i4_to_i8(3*this%natoms),
data, &
1467 ierr, nohead = .
true., fendian = io_binary_is_little_endian())
1468 safe_deallocate_a(data)
1469 iunit = io_open(trim(fullname), this%namespace, action=
'write', position =
'append')
1470 write(iunit,
'(1a)')
''
1473 write(iunit,
'(a,2i9)')
'VERTICES ', this%natoms, 2*this%natoms
1476 do iatom = 1, this%natoms
1477 write(iunit,
'(2i9)') 1, iatom - 1
1480 call io_close(iunit)
1481 safe_allocate(idata(1:2, 1:this%natoms))
1482 do iatom = 1, this%natoms
1484 idata(2, iatom) = iatom - 1
1486 call io_binary_write(io_workpath(fullname, this%namespace), i4_to_i8(2*this%natoms), idata, &
1487 ierr, nohead = .
true., fendian = io_binary_is_little_endian())
1488 safe_deallocate_a(idata)
1489 iunit = io_open(trim(fullname), this%namespace, action=
'write', position =
'append')
1490 write(iunit,
'(1a)')
''
1493 write(iunit,
'(a,i9)')
'POINT_DATA', this%natoms
1494 write(iunit,
'(a)')
'SCALARS element integer'
1495 write(iunit,
'(a)')
'LOOKUP_TABLE default'
1498 do iatom = 1, this%natoms
1499 write(iunit,
'(i9)') nint(this%atom(iatom)%species%get_z())
1502 call io_close(iunit)
1504 safe_allocate(idata(1:this%natoms, 1))
1506 do iatom = 1, this%natoms
1507 idata(iatom, 1) = nint(this%atom(iatom)%species%get_z())
1510 call io_binary_write(io_workpath(fullname, this%namespace), i4_to_i8(this%natoms), idata, &
1511 ierr, nohead = .
true., fendian = io_binary_is_little_endian())
1513 safe_deallocate_a(idata)
1515 iunit = io_open(trim(fullname), this%namespace, action=
'write', position =
'append')
1516 write(iunit,
'(1a)')
''
1519 call io_close(iunit)
1526 class(
ions_t),
intent(in) :: this
1527 real(real64) :: current(this%space%dim)
1534 do iatom = this%atoms_dist%start, this%atoms_dist%end
1535 current = current + p_proton_charge*this%atom(iatom)%species%get_zval()*this%vel(:,iatom)
1538 if (this%atoms_dist%parallel)
then
1539 call comm_allreduce(this%atoms_dist%mpi_grp, current, dim=this%space%dim)
1547 class(
ions_t),
intent(in) :: this
1548 real(real64) :: abs_current(this%space%dim)
1554 abs_current = m_zero
1555 do iatom = this%atoms_dist%start, this%atoms_dist%end
1556 abs_current = abs_current + abs(p_proton_charge*this%atom(iatom)%species%get_zval()*this%vel(:,iatom))
1559 if (this%atoms_dist%parallel)
then
1560 call comm_allreduce(this%atoms_dist%mpi_grp, abs_current, dim=this%space%dim)
1568 type(
ions_t),
intent(inout) :: ions
1574 call distributed_end(ions%atoms_dist)
1576 call ion_interaction_end(ions%ion_interaction)
1578 safe_deallocate_a(ions%atom)
1581 if(
allocated(ions%species))
then
1582 do i = 1, ions%nspecies
1584 if(
associated(ions%species(i)%s))
deallocate(ions%species(i)%s)
1586 deallocate(ions%species)
1590 safe_deallocate_a(ions%pos_displacements)
1591 safe_deallocate_a(ions%vel_displacements)
1593 call charged_particles_end(ions)
1595 safe_deallocate_a(ions%map_symm_atoms)
1596 safe_deallocate_a(ions%inv_map_symm_atoms)
1606 class(
ions_t),
intent(inout) :: ions
1607 type(lattice_vectors_t),
intent(in) :: latt
1608 logical,
intent(in) :: symmetrize
1613 call ions%latt%update(latt%rlattice)
1616 if (symmetrize)
then
1617 call symmetries_update_lattice_vectors(ions%symm, latt, ions%space%dim)
1630 class(
ions_t),
intent(inout) :: ions
1632 integer :: iatom, iop, iatom_symm, dim4symms
1633 real(real64) :: ratom(ions%space%dim)
1637 safe_allocate(ions%map_symm_atoms(1:ions%natoms, 1:ions%symm%nops + ions%symm%nops_nonsymmorphic))
1638 safe_allocate(ions%inv_map_symm_atoms(1:ions%natoms, 1:ions%symm%nops + ions%symm%nops_nonsymmorphic))
1641 dim4symms = min(3, ions%space%dim)
1644 do iop = 1, ions%symm%nops
1645 do iatom = 1, ions%natoms
1647 ratom(1:dim4symms) = symm_op_apply_inv_cart(ions%symm%ops(iop), ions%pos(:, iatom))
1649 ratom(:) = ions%latt%fold_into_cell(ratom(:))
1652 do iatom_symm = 1, ions%natoms
1653 if (all(abs(ratom(:) - ions%pos(:, iatom_symm)) < symprec))
exit
1656 if (iatom_symm > ions%natoms)
then
1657 write(message(1),
'(a,i6)')
'Internal error: could not find symmetric partner for atom number', iatom
1658 write(message(2),
'(a,i3,a)')
'with symmetry operation number ', iop,
'.'
1659 call messages_fatal(2, namespace=ions%namespace)
1662 ions%map_symm_atoms(iatom, iop) = iatom_symm
1663 ions%inv_map_symm_atoms(iatom_symm, iop) = iatom
1668 do iop = 1, ions%symm%nops_nonsymmorphic
1669 do iatom = 1, ions%natoms
1671 ratom(1:dim4symms) = symm_op_apply_inv_cart(ions%symm%non_symmorphic_ops(iop), ions%pos(:, iatom))
1673 ratom(:) = ions%latt%fold_into_cell(ratom(:))
1676 do iatom_symm = 1, ions%natoms
1677 if (all(abs(ratom(:) - ions%pos(:, iatom_symm)) < symprec))
exit
1680 if (iatom_symm > ions%natoms)
then
1681 write(message(1),
'(a,i6)')
'Internal error: could not find symmetric partner for atom number', iatom
1682 write(message(2),
'(a,i3,a)')
'with symmetry operation number ', iop,
'.'
1683 call messages_fatal(2, namespace=ions%namespace)
1686 ions%map_symm_atoms(iatom, ions%symm%nops+iop) = iatom_symm
1687 ions%inv_map_symm_atoms(iatom_symm, ions%symm%nops+iop) = iatom
1702 integer :: iatom, iop, iatom_sym
1703 real(real64) :: ratom(ions%space%dim)
1704 real(real64),
allocatable :: new_pos(:,:)
1708 safe_allocate(new_pos(1:ions%space%dim, 1:ions%natoms))
1710 do iatom = 1, ions%natoms
1711 new_pos(:, iatom) = m_zero
1714 do iop = 1, ions%symm%nops
1715 iatom_sym = ions%inv_map_symm_atoms(iatom, iop)
1716 ratom = symm_op_apply_inv_cart(ions%symm%ops(iop), ions%pos(:, iatom_sym))
1717 ratom = ions%latt%fold_into_cell(ratom)
1718 new_pos(:, iatom) = new_pos(:, iatom) + ratom
1722 do iop = 1, ions%symm%nops_nonsymmorphic
1723 iatom_sym = ions%inv_map_symm_atoms(iatom, iop + ions%symm%nops)
1724 ratom = symm_op_apply_inv_cart(ions%symm%non_symmorphic_ops(iop), ions%pos(:, iatom_sym))
1725 ratom = ions%latt%fold_into_cell(ratom)
1726 new_pos(:, iatom) = new_pos(:, iatom) + ratom
1729 new_pos(:, iatom) = new_pos(:, iatom) / (ions%symm%nops + ions%symm%nops_nonsymmorphic)
1738 class(
ions_t),
intent(inout) :: ions
1740 type(spglibdataset) :: spg_dataset
1741 character(len=11) :: symbol
1742 integer,
allocatable :: site_type(:)
1743 integer :: space_group, ia
1745 if(.not. ions%space%is_periodic())
return
1749 safe_allocate(site_type(1:ions%natoms))
1750 do ia = 1, ions%natoms
1751 site_type(ia) = ions%atom(ia)%species%get_index()
1754 spg_dataset = symmetries_get_spg_dataset(ions%namespace, ions%space, ions%latt, ions%natoms, ions%pos, site_type)
1756 safe_deallocate_a(site_type)
1758 if (spg_dataset%spglib_error /= 0)
then
1763 space_group = spg_dataset%spacegroup_number
1764 symbol = spg_dataset%international_symbol
1766 write(message(1),
'(a, i4)')
'Info: Space group No. ', space_group
1767 write(message(2),
'(2a)')
'Info: International: ', trim(symbol)
1768 call messages_info(2, namespace=ions%namespace)
--------------— axpy ---------------— Constant times a vector plus a vector.
double acos(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
double tanh(double __x) __attribute__((__nothrow__
double atan2(double __y, double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
subroutine rotate(m, angle, dir)
subroutine, public atom_init(this, dim, label, species)
subroutine, public atom_get_species(this, species)
subroutine, public atom_set_species(this, species)
This module contains interfaces for BLAS routines You should not use these routines directly....
This module handles the calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
subroutine, public charged_particles_copy(this, cp_in)
subroutine, public charged_particles_init_interaction_as_partner(partner, interaction)
subroutine, public charged_particles_update_quantity(this, label)
subroutine, public charged_particles_init_interaction(this, interaction)
subroutine, public charged_particles_copy_quantities_to_interaction(partner, interaction)
subroutine, public charged_particles_init(this, np)
The init routine is a module level procedure This has the advantage that different classes can have d...
subroutine, public distributed_nullify(this, total)
subroutine, public distributed_init(this, total, comm, tag, scalapack_compat)
Distribute N instances across M processes of communicator comm
subroutine, public distributed_copy(in, out)
@Brief Create a copy of a distributed instance
real(real64), parameter, public m_huge
real(real64), parameter, public m_zero
character(len= *), parameter, public static_dir
real(real64), parameter, public p_kb
Boltzmann constant in Ha/K.
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
subroutine, public io_mkdir(fname, namespace, parents)
subroutine, public ion_interaction_init_parallelization(this, natoms, mc)
subroutine, public ion_interaction_init(this, namespace, space, natoms)
subroutine ions_init_interaction(this, interaction)
subroutine ions_update_lattice_vectors(ions, latt, symmetrize)
Regenerate the ions information after update of the lattice vectors.
class(ions_t) function, pointer ions_constructor(namespace, print_info, latt_inp)
subroutine ions_fold_atoms_into_cell(this)
subroutine ions_single_mode_displacements(ions, mode, amplitude_pos, amplitude_vel)
apply initial displacements and velocities corresponding to a given phonon mode
real(real64) function, dimension(this%space%dim) ions_global_force(this, time)
subroutine ions_copy(ions_out, ions_in)
real(real64) function, dimension(this%space%dim) ions_current(this)
subroutine ions_update_quantity(this, label)
subroutine ions_generate_mapping_symmetric_atoms(ions)
Given the symmetries of the system, we create a mapping that tell us for each atom and symmetry,...
subroutine ions_write_xyz(this, fname, append, comment, reduce_coordinates)
subroutine ions_init_species(ions, factory, print_info)
subroutine ions_finalize(ions)
subroutine ions_initialize(this)
real(real64) function, dimension(this%space%dim) ions_abs_current(this)
subroutine ions_read_xyz(this, fname, comment)
subroutine ions_write_vtk_geometry(this, filename, ascii)
subroutine ions_rotate(this, from, from2, to)
subroutine ions_translate(this, xx)
subroutine ions_symmetrize_atomic_coord(ions)
Symmetrizes atomic coordinates by applying all symmetries.
subroutine ions_print_spacegroup(ions)
Prints the spacegroup of the system for periodic systems.
real(real64) function ions_min_distance(this, real_atoms_only)
subroutine ions_init_random_displacements(ions, T)
create random displacements for positions and velocities
real(real64) function, dimension(this%space%dim) ions_dipole(this, mask)
subroutine ions_write_bild_forces_file(this, dir, fname)
subroutine ions_write_crystal(this, dir)
This subroutine creates a crystal by replicating the geometry and writes the result to dir.
logical function ions_has_time_dependent_species(this)
subroutine ions_partition(this, mc)
subroutine ions_init_interaction_as_partner(partner, interaction)
subroutine ions_copy_quantities_to_interaction(partner, interaction)
subroutine ions_write_poscar(this, fname, comment)
Writes the positions of the ions in POSCAR format.
real(real64) function ions_val_charge(this, mask)
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
character(len=512), private msg
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)
This module handles the communicators for the various parallelization strategies.
logical function, public parse_is_defined(namespace, name)
This module provides a class for (classical) phonon modes.
integer, parameter, public read_coords_reduced
subroutine, public read_coords_init(gf)
integer, parameter, public xyz_flags_move
subroutine, public read_coords_end(gf)
subroutine, public read_coords_read(what, gf, space, namespace)
subroutine, public species_factory_init(factory, namespace)
subroutine, public species_factory_end(factory)
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_t), public unit_amu
Mass in atomic mass units (AKA Dalton).
type(unit_system_t), public units_out
integer, parameter, public wigner_q
integer, parameter, public wigner_p
abstract interaction class
surrogate interaction class to avoid circular dependencies between modules.
Stores all communicators and groups.
This class describes phonon modes, which are specified by their frequencies and eigenvectors.
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...
Class describing a Wigner distribution for sampling initial conditions in multi-trajectory Ehrenfest ...