Octopus
ions.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module ions_oct_m
23 use atom_oct_m
24 use blas_oct_m
27 use comm_oct_m
28 use iso_c_binding
29 use debug_oct_m
31 use global_oct_m
35 use io_oct_m
37 use, intrinsic :: iso_fortran_env
41 use math_oct_m
44 use mpi_oct_m
46 use parser_oct_m
51 use space_oct_m
54 use spglib_f08
58 use unit_oct_m
62
63 implicit none
64
65 private
66 public :: ions_t
67
68 type, extends(charged_particles_t) :: ions_t
69 ! Components are public by default
70
71 type(lattice_vectors_t) :: latt
72
73 integer :: natoms
74 type(atom_t), allocatable :: atom(:)
75
76 type(symmetries_t) :: symm
77
78 type(distributed_t) :: atoms_dist
79
80 integer, allocatable :: map_symm_atoms(:,:)
81 integer, allocatable :: inv_map_symm_atoms(:,:)
82
83 real(real64), allocatable :: equilibrium_pos(:,:)
84 ! !! for multitrajectory runs.
85 ! !!
86 ! !! This is necessary, as the displacements are added before we
87 ! !! generate the grid, but the grid should always be constructed
88 ! !! for the equilibrium positions.
89 real(real64), allocatable :: pos_displacements(:,:)
90 real(real64), allocatable :: vel_displacements(:,:)
91
92 ! Information about the species
93 integer :: nspecies
94 class(species_wrapper_t), allocatable :: species(:)
95 logical :: only_user_def
96 logical, private :: species_time_dependent
97
98 logical :: force_total_enforce
99 type(ion_interaction_t) :: ion_interaction
100
102 logical, private :: apply_global_force
103 type(tdf_t), private :: global_force_function
104 contains
105 procedure :: copy => ions_copy
106 generic :: assignment(=) => copy
107 procedure :: partition => ions_partition
108 procedure :: init_interaction => ions_init_interaction
109 procedure :: initialize => ions_initialize
110 procedure :: update_quantity => ions_update_quantity
111 procedure :: init_interaction_as_partner => ions_init_interaction_as_partner
112 procedure :: copy_quantities_to_interaction => ions_copy_quantities_to_interaction
113 procedure :: fold_atoms_into_cell => ions_fold_atoms_into_cell
114 procedure :: min_distance => ions_min_distance
115 procedure :: has_time_dependent_species => ions_has_time_dependent_species
116 procedure :: val_charge => ions_val_charge
117 procedure :: dipole => ions_dipole
118 procedure :: translate => ions_translate
119 procedure :: rotate => ions_rotate
120 procedure :: global_force => ions_global_force
121 procedure :: write_xyz => ions_write_xyz
122 procedure :: read_xyz => ions_read_xyz
123 procedure :: write_poscar => ions_write_poscar
124 procedure :: write_crystal => ions_write_crystal
125 procedure :: write_bild_forces_file => ions_write_bild_forces_file
126 procedure :: write_vtk_geometry => ions_write_vtk_geometry
127 procedure :: update_lattice_vectors => ions_update_lattice_vectors
128 procedure :: symmetrize_atomic_coord => ions_symmetrize_atomic_coord
129 procedure :: generate_mapping_symmetric_atoms => ions_generate_mapping_symmetric_atoms
130 procedure :: print_spacegroup => ions_print_spacegroup
131 procedure :: current => ions_current
132 procedure :: abs_current => ions_abs_current
133 procedure :: init_random_displacements => ions_init_random_displacements
134 procedure :: single_mode_displacements => ions_single_mode_displacements
135 final :: ions_finalize
136 end type ions_t
137
138 interface ions_t
139 procedure ions_constructor
140 end interface ions_t
141
142contains
143
144 ! ---------------------------------------------------------
145 function ions_constructor(namespace, print_info, latt_inp) result(ions)
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
150
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
160 real(real64) :: T
161 integer :: displacement_mode
162 real(real64) :: amp_pos, amp_vel
164 logical :: in_ensemble
165
166 push_sub_with_profile(ions_constructor)
167
168 allocate(ions)
170 ions%namespace = namespace
172 ions%space = space_t(namespace)
173 in_ensemble = parse_is_defined(namespace, "NumberOfReplicas")
174
176 call species_factory_init(factory, namespace)
177
178 ! initialize geometry
179 call read_coords_init(xyz)
180
181 ! load positions of the atoms
182 call read_coords_read('Coordinates', xyz, ions%space, namespace)
183
184 if (xyz%n < 1) then
185 message(1) = "Coordinates have not been defined."
186 call messages_fatal(1, namespace=namespace)
187 end if
189 ! Initialize parent class
190 call charged_particles_init(ions, xyz%n)
192 ! copy information from xyz to ions
193 ions%natoms = xyz%n
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)
198 if (bitand(xyz%flags, xyz_flags_move) /= 0) then
199 ions%fixed(ia) = .not. xyz%atom(ia)%move
200 end if
201 end do
203 if (allocated(xyz%latvec)) then
204 ! Build lattice vectors from the XSF input
205 ions%latt = lattice_vectors_t(namespace, ions%space, xyz%latvec)
206 else
207 ! Build lattice vectors from input file
208 ions%latt = lattice_vectors_t(namespace, ions%space)
209 end if
211 ! Convert coordinates to Cartesian in case we have reduced coordinates
212 if (xyz%source == read_coords_reduced) then
213 do ia = 1, ions%natoms
214 ions%pos(:, ia) = ions%latt%red_to_cart(ions%pos(:, ia))
215 end do
216 end if
220 ! Save the positions read from the input before applying deviations from the initial_input
221 safe_allocate_source_a(ions%equilibrium_pos, ions%pos)
223 call ions_init_species(ions, factory, print_info=print_info)
225 ! Set the masses and charges. This needs to be done after initializing the species.
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()
229 end do
231 !%Variable InitialDisplacementMode
232 !%Type integer
233 !%Default 0
234 !%Section System
235 !%Description
236 !% When this variable is set and non-zero, the phonon modes file will be read, and ionic
237 !% positions and velocities will be displaced according to a single phonon mode.
238 !% The amplitudes can be defined with the variables:
239 !% - InitialDisplacementsAmplitudePos
240 !% - InitialDisplacementsAmplitudeVel
241 !%End
242 call parse_variable(namespace, 'InitialDisplacementMode', 0, displacement_mode)
243
244 !%Variable InitialDisplacementAmplitudePos
245 !%Type float
246 !%Default 1.0
247 !%Section System
248 !%Description
249 !% Amplitude for initial position displacement according to the phonon mode, defined by InitialDisplacementMode.
250 !% The value corresponds to the normal distributed canonical coordinates.
251 !%End
252 call parse_variable(namespace, 'InitialDisplacementAmplitudePos', 1.0_real64, amp_pos)
253
254 !%Variable InitialDisplacementAmplitudeVel
255 !%Type float
256 !%Default 1.0
257 !%Section System
258 !%Description
259 !% Amplitude for initial velocity displacement according to the phonon mode, defined by InitialDisplacementMode.
260 !% The value corresponds to the normal distributed canonical coordinates.
261 !%End
262 call parse_variable(namespace, 'InitialDisplacementAmplitudeVel', 1.0_real64, amp_vel)
263
264 if (displacement_mode>0) then
265
266 message(1) = "Random displacements are disabled when InitialDisplacementMode > 0."
267 call messages_warning(1, namespace=namespace)
268
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))
271
272 call ions%single_mode_displacements(displacement_mode, amp_pos, amp_vel)
273
274 ! apply initial displacements
275 ions%pos = ions%pos + ions%pos_displacements
276
277 ! note that the velocity displacement needs to be applied in ions%initialize()
278
279 end if
280
281 ! Now we can apply displacements to the positions, if requested.
282 ! The random displacements are disabled if a single mode is specified,
283 if(in_ensemble .and. displacement_mode==0) then
284
285 !%Variable EnsembleTemperature
286 !%Type float
287 !%Default 0
288 !%Section Multi-Trajectory
289 !%Description
290 !% The temperature for an ensemble calculation in Kelvin
291 !%End
292 call parse_variable(namespace, 'EnsembleTemperature', 0.0_real64, t)
293
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)
297
298 ! apply initial displacements
299 ions%pos = ions%pos + ions%pos_displacements
300
301 ! note that the velocity displacement needs to be applied in ions%initialize()
302
303 end if
304
306 call distributed_nullify(ions%atoms_dist, ions%natoms)
307
308 call messages_obsolete_variable(namespace, 'PDBClassical')
309
310 if (present(latt_inp)) then
311 ! The lattice as read from the input might be needed by some other part of the code, so we save it
312 latt_inp = ions%latt
313 end if
314
315 ! Now that we have processed the atomic coordinates, we renormalize the
316 ! lattice parameters along the non-periodic dimensions
317 if (ions%space%has_mixed_periodicity()) then
318 safe_allocate(factor(ions%space%dim))
319 do idir = 1, ions%space%periodic_dim
320 factor(idir) = m_one
321 end do
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))
324 end do
325 call ions%latt%scale(factor)
326 safe_deallocate_a(factor)
327 end if
328
329 ! Check that atoms are not too close
330 if (ions%natoms > 1) then
331 mindist = ions_min_distance(ions, real_atoms_only = .false.)
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 = ", &
336 units_from_atomic(units_out%length, mindist), trim(units_abbrev(units_out%length))
337 call messages_warning(3, namespace=namespace)
338
339 ! then write out the geometry, whether asked for or not in Output variable
340 call io_mkdir(static_dir, namespace)
341 call ions%write_xyz(trim(static_dir)//'/geometry')
342 end if
343
344 if (ions_min_distance(ions, real_atoms_only = .true.) < threshold) then
345 message(1) = "It cannot be correct to run with physical atoms so close."
346 call messages_fatal(1, namespace=namespace)
347 end if
348 end if
349
350 !Initialize symmetries
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)
355 type is(jellium_slab_t)
356 spherical_site(ia) = .false.
358 spherical_site(ia) = .false.
359 type is(species_from_file_t)
360 spherical_site(ia) = .false.
362 spherical_site(ia) = .false.
363 class default
364 spherical_site(ia) = .true.
365 end select
366
367 site_type(ia) = ions%atom(ia)%species%get_index()
368 end do
369
370 ions%symm = symmetries_t(ions%namespace, ions%space, ions%latt, ions%natoms, ions%pos, site_type, spherical_site)
371
372 safe_deallocate_a(spherical_site)
373 safe_deallocate_a(site_type)
374
375 ! Generate the mapping of symmetric atoms
376 call ions%generate_mapping_symmetric_atoms()
377
378 call ion_interaction_init(ions%ion_interaction, namespace, ions%space, ions%natoms)
379
380 !%Variable ForceTotalEnforce
381 !%Type logical
382 !%Default no
383 !%Section Hamiltonian
384 !%Description
385 !% (Experimental) If this variable is set to "yes", then the sum
386 !% of the total forces will be enforced to be zero.
387 !%End
388 call parse_variable(namespace, 'ForceTotalEnforce', .false., ions%force_total_enforce)
389 if (ions%force_total_enforce) call messages_experimental('ForceTotalEnforce', namespace=namespace)
390
391 !%Variable TDGlobalForce
392 !%Type string
393 !%Section Time-Dependent
394 !%Description
395 !% If this variable is set, a global time-dependent force will be
396 !% applied to the ions in the x direction during a time-dependent
397 !% run. This variable defines the base name of the force, that
398 !% should be defined in the <tt>TDFunctions</tt> block. This force
399 !% does not affect the electrons directly.
400 !%End
401
402 if (parse_is_defined(namespace, 'TDGlobalForce')) then
403
404 ions%apply_global_force = .true.
405
406 call parse_variable(namespace, 'TDGlobalForce', 'none', function_name)
407 call tdf_read(ions%global_force_function, namespace, trim(function_name), ierr)
408
409 if (ierr /= 0) then
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.")
412 call messages_fatal(namespace=namespace)
413 end if
414
415 else
416
417 ions%apply_global_force = .false.
418
419 end if
420
421 call species_factory_end(factory)
422
423 pop_sub_with_profile(ions_constructor)
424 end function ions_constructor
425
426 ! ---------------------------------------------------------
427 subroutine ions_init_species(ions, factory, print_info)
428 type(ions_t), intent(inout) :: ions
429 type(species_factory_t), intent(in) :: factory
430 logical, optional, intent(in) :: print_info
431
432 logical :: print_info_, spec_user_defined
433 integer :: i, j, k, ispin
434 class(species_t), pointer :: spec
435
436 push_sub_with_profile(ions_init_species)
437
438 print_info_ = .true.
439 if (present(print_info)) then
440 print_info_ = print_info
441 end if
442 ! First, count the species
443 ions%nspecies = 0
444 atoms1: do i = 1, ions%natoms
445 do j = 1, i - 1
446 if (atom_same_species(ions%atom(j), ions%atom(i))) cycle atoms1
447 end do
448 ions%nspecies = ions%nspecies + 1
449 end do atoms1
450
451 ! Allocate the species structure.
452 allocate(ions%species(1:ions%nspecies))
453
454 ! Now, read the data.
455 k = 0
456 ions%only_user_def = .true.
457 atoms2: do i = 1, ions%natoms
458 do j = 1, i - 1
459 if (atom_same_species(ions%atom(j), ions%atom(i))) cycle atoms2
460 end do
461 k = k + 1
462 ions%species(k)%s => factory%create_from_input(ions%namespace, ions%atom(j)%get_label(), k)
463 ! these are the species which do not represent atoms
464 select type(spec => ions%species(k)%s)
465 class is(jellium_t)
466
467 class default
468 ions%only_user_def = .false.
469 end select
470
471 select type(spec => ions%species(k)%s)
472 type is(pseudopotential_t)
473 if (ions%space%dim /= 3) then
474 message(1) = "Pseudopotentials may only be used with Dimensions = 3."
475 call messages_fatal(1, namespace=ions%namespace)
476 end if
477
478 type is(jellium_slab_t)
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"
481 call messages_fatal(1, namespace=ions%namespace)
482 end if
483 end select
484
485 end do atoms2
486
487 ! Reads the spin components. This is read here, as well as in states_init,
488 ! to be able to pass it to the pseudopotential initializations subroutine.
489 call parse_variable(ions%namespace, 'SpinComponents', 1, ispin)
490 if (.not. varinfo_valid_option('SpinComponents', ispin)) call messages_input_error(ions%namespace, 'SpinComponents')
491 ispin = min(2, ispin)
492
493 if (print_info_) then
494 call messages_print_with_emphasis(msg="Species", namespace=ions%namespace)
495 end if
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_)
499 end do
500 if (print_info_) then
501 call messages_print_with_emphasis(namespace=ions%namespace)
502 end if
503
504 !%Variable SpeciesTimeDependent
505 !%Type logical
506 !%Default no
507 !%Section System::Species
508 !%Description
509 !% When this variable is set, the potential defined in the block <tt>Species</tt> is calculated
510 !% and applied to the Hamiltonian at each time step. You must have at least one <tt>species_user_defined</tt>
511 !% type of species to use this.
512 !%End
513 call parse_variable(ions%namespace, 'SpeciesTimeDependent', .false., ions%species_time_dependent)
514 ! we must have at least one user defined species in order to have time dependency
515 do i = 1,ions%nspecies
516 select type(spec=>ions%species(i)%s)
518 spec_user_defined = .true.
519 end select
520 end do
521 if (ions%species_time_dependent .and. .not. spec_user_defined) then
522 call messages_input_error(ions%namespace, 'SpeciesTimeDependent')
523 end if
524
525 ! assign species
526 do i = 1, ions%natoms
527 do j = 1, ions%nspecies
528 if (atom_same_species(ions%atom(i), ions%species(j)%s)) then
529 call atom_set_species(ions%atom(i), ions%species(j)%s)
530 exit
531 end if
532 end do
533 end do
534
535 pop_sub_with_profile(ions_init_species)
536 end subroutine ions_init_species
537
538 !--------------------------------------------------------------
539 subroutine ions_copy(ions_out, ions_in)
540 class(ions_t), intent(out) :: ions_out
541 class(ions_t), intent(in) :: ions_in
542
543 push_sub(ions_copy)
544
545 call charged_particles_copy(ions_out, ions_in)
546
547 ions_out%latt = ions_in%latt
548
549 ions_out%natoms = ions_in%natoms
550 safe_allocate(ions_out%atom(1:ions_out%natoms))
551 ions_out%atom = ions_in%atom
552
553 ions_out%nspecies = ions_in%nspecies
554 allocate(ions_out%species(1:ions_out%nspecies))
555 ions_out%species = ions_in%species
556
557 ions_out%only_user_def = ions_in%only_user_def
558
559 call distributed_copy(ions_in%atoms_dist, ions_out%atoms_dist)
560
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
565
566
567 pop_sub(ions_copy)
568 end subroutine ions_copy
569
570 ! ---------------------------------------------------------
571 subroutine ions_partition(this, mc)
572 class(ions_t), intent(inout) :: this
573 type(multicomm_t), intent(in) :: mc
574
575 push_sub(ions_partition)
576
577 call distributed_init(this%atoms_dist, this%natoms, mc%group_comm(p_strategy_states), "atoms")
578
579 call ion_interaction_init_parallelization(this%ion_interaction, this%natoms, mc)
580
581 pop_sub(ions_partition)
582 end subroutine ions_partition
583
584 ! ---------------------------------------------------------
585 subroutine ions_init_interaction(this, interaction)
586 class(ions_t), target, intent(inout) :: this
587 class(interaction_t), intent(inout) :: interaction
588
589 push_sub(ions_init_interaction)
590
591 select type (interaction)
592 class default
593 call charged_particles_init_interaction(this, interaction)
594 end select
595
596 pop_sub(ions_init_interaction)
597 end subroutine ions_init_interaction
598
599 ! ---------------------------------------------------------
607 !
608 subroutine ions_init_random_displacements(ions, T)
609 class(ions_t), intent(inout) :: ions
610 real(real64), intent(in) :: T
611
612
613 integer(int64) :: seed
614
615 type(phonon_modes_t) :: phonons
616 type(wigner_distribution_t) :: wigner
617
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(:)
624
625
626 real(real64) :: beta_half
627 real(real64), allocatable :: pos_displacements(:), vel_displacements(:)
628 real(real64), allocatable :: l_m(:), prefactor(:)
629
630 integer :: imode, idim, iatom, i, iline
631 integer :: num_real_modes
632
633 real(real64), parameter :: low_temperature_tolerance = 1.0e-6_real64
635
637
638 ! create a unique seed for each replica, based on a hash of the full (unique) namespace string.
639 seed = ions%namespace%get_hash32()
640
641 message(1) = "Create initial random displacements for ions."
642 call messages_info(1, namespace = ions%namespace)
643
644 write(message(1), '("namespace = ",A)') ions%namespace%get()
645 write(message(2), '("seed = ",I0)') seed
646 call messages_info(2, namespace = ions%namespace, debug_only=.true.)
647
648
649 ! initialize the phonons (read info from file)
650 call phonons%init(ions%namespace, ions%space%dim, ions%natoms, ions%space%periodic_dim > 0)
651
652 num_real_modes = phonons%num_modes
653
654 if (num_real_modes == 0) then
655
656 ions%pos_displacements = m_zero
657 ions%vel_displacements = m_zero
658
660 return
661 end if
662
663
664 ! initialize wigner distribution for phonons%num_modes modes and with given seed
665 call wigner%init(num_real_modes, seed)
667 ! set up widths and shifts for the Wigner function
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))
674
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)))
678
679 safe_allocate(l_m(1:num_real_modes))
681
682 if (t < low_temperature_tolerance) then
683 ! TODO: Implement proper low T limit
684 sigmap = 1.0_real64/2.0_real64
685 sigmaq = 1.0_real64/2.0_real64
686 else
687 beta_half = m_one / (2 * p_kb * t)
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)))
690 end if
691 mup = m_zero
692 muq = m_zero
693
694 l_m = sqrt(2.0_real64/(unit_amu%factor * phonons%frequencies)) * phonons%alpha
695
696 i = 1
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)
700 i=i+1
701 end do
702 end do
704 ! get generalized mode coords from Wigner function
705 genq = wigner%get(sigmaq, muq, wigner_q) * l_m
706 genp = wigner%get(sigmap, mup, wigner_p) / (l_m * unit_amu%factor) ! TODO: Check!!
707
708 ! construct initial displacements and velocities for ions (velocities must be applied later)
709
710 ions%pos_displacements = m_zero
711 ions%vel_displacements = m_zero
712
713 ! implement eq.(9) from Kevins paper for periodic systems
714
715 do imode = 1, num_real_modes
716
717 pos_displacements = prefactor(:) * phonons%eigenvectors(:, imode) * genq(imode)
718 vel_displacements = prefactor(:) * phonons%eigenvectors(:, imode) * genp(imode)
719
720 iline = 1
721 write(message(iline), '("Displacements for mode ",I4)') imode
722 do iatom=1, ions%natoms
723 iline = iline+1
724 write(message(iline), '(3E15.5)') pos_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
725 end do
726 iline = iline+1
727 write(message(iline), '("Velocities for mode ",I4)') imode
728 do iatom=1, ions%natoms
729 iline = iline+1
730 write(message(iline), '(3E15.5)') vel_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
731 end do
732
733 call messages_info(iline, namespace=ions%namespace)
734
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)
737
738 end do
739
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)
747
748 safe_deallocate_a(prefactor)
749 safe_deallocate_a(pos_displacements)
750 safe_deallocate_a(vel_displacements)
751
753
754 end subroutine ions_init_random_displacements
755
756
760 !
761 subroutine ions_single_mode_displacements(ions, mode, amplitude_pos, amplitude_vel)
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
766
767 type(phonon_modes_t) :: phonons
768 integer :: num_real_modes, i, iatom, idim
769 real(real64) :: l_m, genP, genQ
770
771 real(real64), allocatable :: pos_displacements(:), vel_displacements(:), prefactor(:)
772
773
775
776 write(message(1), '("Create displacements for mode ", I4)') mode
777 call messages_info(1, namespace=ions%namespace)
778
779 ! initialize the phonons (read info from file)
780 call phonons%init(ions%namespace, ions%space%dim, ions%natoms, ions%space%periodic_dim > 0)
781
782 num_real_modes = phonons%num_modes
783
784 if (mode > num_real_modes) then
785 write(message(1), '("Requested mode ",I0," exceeds number of available modes (",I0,")")') mode, num_real_modes
786 call messages_fatal(1, namespace=ions%namespace)
787
789 return
790 end if
791
792 ions%pos_displacements = m_zero
793 ions%vel_displacements = m_zero
794
795 l_m = sqrt(2.0_real64/(unit_amu%factor * phonons%frequencies(mode))) * phonons%alpha(mode)
796
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)))
800
801 i = 1
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)
805 i=i+1
806 end do
807 end do
808
809 ! get generalized mode coords from Wigner function
810 genq = amplitude_pos * l_m
811 genp = amplitude_vel / (l_m * unit_amu%factor)! TODO: Check!!
812
813 pos_displacements = prefactor(:) * phonons%eigenvectors(:, mode) * genq
814 vel_displacements = prefactor(:) * phonons%eigenvectors(:, mode) * genp
815
816 write(message(1), '("Displacements for mode ",I4)') mode
817 call messages_info(1, namespace=ions%namespace)
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)
820 call messages_info(1, namespace=ions%namespace)
821 end do
822 write(message(1), '("Velocities for mode ",I4)') mode
823 call messages_info(1, namespace=ions%namespace)
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)
826 call messages_info(1, namespace=ions%namespace)
827 end do
828
829
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)
832
833 safe_deallocate_a(pos_displacements)
834 safe_deallocate_a(vel_displacements)
835 safe_deallocate_a(prefactor)
836
838
839 end subroutine ions_single_mode_displacements
840
841 ! ---------------------------------------------------------
842 subroutine ions_initialize(this)
843 class(ions_t), intent(inout) :: this
844
845 push_sub(ions_initialize)
846
847 ! At this point, we need to apply initial velocities from the random displacements.
848 if(allocated(this%vel_displacements)) then
849 this%vel = this%vel + this%vel_displacements
850 end if
851
852 pop_sub(ions_initialize)
853 end subroutine ions_initialize
854
855 ! ---------------------------------------------------------
856 subroutine ions_update_quantity(this, label)
857 class(ions_t), intent(inout) :: this
858 character(len=*), intent(in) :: label
859
860 push_sub(ions_update_quantity)
861
862 select case (label)
863 case default
864 ! Other quantities should be handled by the parent class
865 call charged_particles_update_quantity(this, label)
866 end select
867
868 pop_sub(ions_update_quantity)
869 end subroutine ions_update_quantity
870
871 ! ---------------------------------------------------------
872 subroutine ions_init_interaction_as_partner(partner, interaction)
873 class(ions_t), intent(in) :: partner
874 class(interaction_surrogate_t), intent(inout) :: interaction
875
877
878 select type (interaction)
879 class default
880 call charged_particles_init_interaction_as_partner(partner, interaction)
881 end select
882
885
886 ! ---------------------------------------------------------
887 subroutine ions_copy_quantities_to_interaction(partner, interaction)
888 class(ions_t), intent(inout) :: partner
889 class(interaction_surrogate_t), intent(inout) :: interaction
890
892
893 select type (interaction)
894 class default
895 ! Other interactions should be handled by the parent class
896 call charged_particles_copy_quantities_to_interaction(partner, interaction)
897 end select
898
901
902 ! ---------------------------------------------------------
903 subroutine ions_fold_atoms_into_cell(this)
904 class(ions_t), intent(inout) :: this
905
906 integer :: iatom
907
909
910 do iatom = 1, this%natoms
911 this%pos(:, iatom) = this%latt%fold_into_cell(this%pos(:, iatom))
912 end do
913
915 end subroutine ions_fold_atoms_into_cell
916
917 ! ---------------------------------------------------------
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
921
922 integer :: iatom, jatom, idir
923 real(real64) :: xx(this%space%dim)
924 logical :: real_atoms_only_
925 class(species_t), pointer :: species
926
927 if (this%natoms == 1 .and. .not. this%space%is_periodic()) then
928 rmin = m_huge
929 return
930 end if
931
932 push_sub(ions_min_distance)
933
934 real_atoms_only_ = optional_default(real_atoms_only, .false.)
935
936 ! Without this line, valgrind complains about a conditional jump on uninitialized variable,
937 ! as atom_get_species has an intent(out) that causes a call to the finalizer (with Ifort)
938 nullify(species)
939
940 rmin = huge(rmin)
941 do iatom = 1, this%natoms
942 call atom_get_species(this%atom(iatom), species)
943 select type(species)
944 class is(jellium_t)
945 if (real_atoms_only_) cycle
946 end select
947 do jatom = iatom + 1, this%natoms
948 call atom_get_species(this%atom(jatom), species)
949 select type(species)
950 class is(jellium_t)
951 if (real_atoms_only_) cycle
952 end select
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
957 xx(idir) = xx(idir) - floor(xx(idir) + m_half)
958 end do
959 xx = this%latt%red_to_cart(xx)
960 end if
961 rmin = min(norm2(xx), rmin)
962 end do
963 end do
964
965 if (.not. (this%only_user_def .and. real_atoms_only_)) then
966 ! what if the nearest neighbors are periodic images?
967 do idir = 1, this%space%periodic_dim
968 rmin = min(rmin, norm2(this%latt%rlattice(:,idir)))
969 end do
970 end if
971
972 ! To avoid numerical instabilities, we round this to 6 digits only
973 if(rmin < huge(rmin)/1e6_real64) then
974 rmin = anint(rmin*1e6_real64)*1.0e-6_real64
975 end if
976
977 pop_sub(ions_min_distance)
978 end function ions_min_distance
979
980 ! ---------------------------------------------------------
981 logical function ions_has_time_dependent_species(this) result(time_dependent)
982 class(ions_t), intent(in) :: this
983
985
986 time_dependent = this%species_time_dependent
987
990
991 ! ---------------------------------------------------------
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(:)
995
996 push_sub(ions_val_charge)
997
998 if (present(mask)) then
999 val_charge = -sum(this%charge, mask=mask)
1000 else
1001 val_charge = -sum(this%charge)
1002 end if
1003
1004 pop_sub(ions_val_charge)
1005 end function ions_val_charge
1006
1007 ! ---------------------------------------------------------
1008 function ions_dipole(this, mask) result(dipole)
1009 class(ions_t), intent(in) :: this
1010 logical, optional, intent(in) :: mask(:)
1011 real(real64) :: dipole(this%space%dim)
1012
1013 integer :: ia
1014
1015 push_sub(ions_dipole)
1016
1017 dipole = m_zero
1018 do ia = 1, this%natoms
1019 if (present(mask)) then
1020 if (.not. mask(ia)) cycle
1021 end if
1022 dipole = dipole + this%charge(ia)*this%pos(:, ia)
1023 end do
1024 dipole = p_proton_charge*dipole
1025
1026 pop_sub(ions_dipole)
1027 end function ions_dipole
1028
1029 ! ---------------------------------------------------------
1030 subroutine ions_translate(this, xx)
1031 class(ions_t), intent(inout) :: this
1032 real(real64), intent(in) :: xx(this%space%dim)
1033
1034 integer :: iatom
1035
1036 push_sub(ions_translate)
1037
1038 do iatom = 1, this%natoms
1039 this%pos(:, iatom) = this%pos(:, iatom) - xx
1040 end do
1041
1042 pop_sub(ions_translate)
1043 end subroutine ions_translate
1044
1045 ! ---------------------------------------------------------
1046 subroutine ions_rotate(this, from, from2, to)
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)
1051
1052 integer :: iatom
1053 real(real64) :: m1(3, 3), m2(3, 3)
1054 real(real64) :: m3(3, 3), f2(3), per(3)
1055 real(real64) :: alpha, r
1056
1057 push_sub(ions_rotate)
1058
1059 if (this%space%dim /= 3) then
1060 call messages_not_implemented("ions_rotate in other than 3 dimensions", namespace=this%namespace)
1061 end if
1062
1063 ! initialize matrices
1064 m1 = diagonal_matrix(3, m_one)
1065
1066 ! rotate the to-axis to the z-axis
1067 if (abs(to(2)) > 1d-150) then
1068 alpha = atan2(to(2), to(1))
1069 call rotate(m1, alpha, 3)
1070 end if
1071 alpha = atan2(norm2(to(1:2)), to(3))
1072 call rotate(m1, -alpha, 2)
1073
1074 ! get perpendicular to z and from
1075 f2 = matmul(m1, from)
1076 per(1) = -f2(2)
1077 per(2) = f2(1)
1078 per(3) = m_zero
1079 r = norm2(per)
1080 if (r > m_zero) then
1081 per = per/r
1082 else
1083 per(2) = m_one
1084 end if
1085
1086 ! rotate perpendicular axis to the y-axis
1087 m2 = diagonal_matrix(3, m_one)
1088 alpha = atan2(per(1), per(2))
1089 call rotate(m2, -alpha, 3)
1090
1091 ! rotate from => to (around the y-axis)
1092 m3 = diagonal_matrix(3, m_one)
1093 alpha = acos(min(sum(from*to), m_one))
1094 call rotate(m3, -alpha, 2)
1095
1096 ! join matrices
1097 m2 = matmul(transpose(m2), matmul(m3, m2))
1098
1099 ! rotate around the z-axis to get the second axis
1100 per = matmul(m2, matmul(m1, from2))
1101 alpha = atan2(per(1), per(2))
1102 call rotate(m2, -alpha, 3) ! second axis is now y
1104 ! get combined transformation
1105 m1 = matmul(transpose(m1), matmul(m2, m1))
1106
1107 ! now transform the coordinates
1108 ! it is written in this way to avoid what I consider a bug in the Intel compiler
1109 do iatom = 1, this%natoms
1110 f2 = this%pos(:, iatom)
1111 this%pos(:, iatom) = matmul(m1, f2)
1112 end do
1113
1114 pop_sub(ions_rotate)
1115 contains
1116
1117 ! ---------------------------------------------------------
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
1122
1123 real(real64) :: aux(3, 3), ca, sa
1124
1126
1127 ca = cos(angle)
1128 sa = sin(angle)
1129
1130 aux = m_zero
1131 select case (dir)
1132 case (1)
1133 aux(1, 1) = m_one
1134 aux(2, 2) = ca
1135 aux(3, 3) = ca
1136 aux(2, 3) = sa
1137 aux(3, 2) = -sa
1138 case (2)
1139 aux(2, 2) = m_one
1140 aux(1, 1) = ca
1141 aux(3, 3) = ca
1142 aux(1, 3) = sa
1143 aux(3, 1) = -sa
1144 case (3)
1145 aux(3, 3) = m_one
1146 aux(1, 1) = ca
1147 aux(2, 2) = ca
1148 aux(1, 2) = sa
1149 aux(2, 1) = -sa
1150 end select
1151
1152 m = matmul(aux, m)
1153
1154 pop_sub(ions_rotate.rotate)
1155 end subroutine rotate
1156
1157 end subroutine ions_rotate
1158
1159 ! ---------------------------------------------------------
1160 function ions_global_force(this, time) result(force)
1161 class(ions_t), intent(in) :: this
1162 real(real64), intent(in) :: time
1163 real(real64) :: force(this%space%dim)
1164
1165 push_sub(ions_global_force)
1166
1167 force = m_zero
1168
1169 if (this%apply_global_force) then
1170 force(1) = tdf(this%global_force_function, time)
1171 end if
1172
1173 pop_sub(ions_global_force)
1174 end function ions_global_force
1175
1176 ! ---------------------------------------------------------
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
1183
1184 integer :: iatom, idim, iunit
1185 character(len=6) position
1186 character(len=19) :: frmt
1187 real(real64) :: red(this%space%dim)
1188
1189 if (.not. mpi_world%is_root()) return
1190
1191 push_sub(ions_write_xyz)
1192
1193 position = 'asis'
1194 if (present(append)) then
1195 if (append) position = 'append'
1196 end if
1197 if(.not.optional_default(reduce_coordinates, .false.)) then
1198 iunit = io_open(trim(fname)//'.xyz', this%namespace, action='write', position=position)
1199 else
1200 iunit = io_open(trim(fname)//'.xyz_red', this%namespace, action='write', position=position)
1201 end if
1202
1203 write(iunit, '(i4)') this%natoms
1204 if (present(comment)) then
1205 write(iunit, '(1x,a)') comment
1206 else
1207 write(iunit, '(1x,a,a)') 'units: ', trim(units_abbrev(units_out%length_xyz_file))
1208 end if
1209
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)
1215 else
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)
1218 end if
1219 end do
1220 call io_close(iunit)
1221
1222 pop_sub(ions_write_xyz)
1223 end subroutine ions_write_xyz
1224
1229 !
1230 subroutine ions_write_poscar(this, fname, comment)
1231 class(ions_t), intent(in) :: this
1232 character(len=*), optional, intent(in) :: fname
1233 character(len=*), optional, intent(in) :: comment
1234
1235 integer :: iatom, idim, iunit
1236 character(len=:), allocatable :: fname_, comment_
1237 character(len=10) :: format_string
1238
1239 push_sub(ions_write_poscar)
1240
1241 if (.not. mpi_world%is_root()) then
1242 pop_sub(ions_write_poscar)
1243 return
1244 end if
1245
1246 comment_ = optional_default(comment, "")
1247 fname_ = optional_default(fname, "POSCAR")
1248
1249 iunit = io_open(trim(fname_), this%namespace, action='write')
1250
1251 write(iunit, '(A)') comment_ ! mandatory comment
1252 write(iunit, '("1.0")') ! scaling factor: we use 1 as lattice vectors are absolute
1253
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))
1257 end do
1258
1259 do iatom = 1, this%natoms
1260 write(iunit, '(A)', advance='NO') trim(this%atom(iatom)%label)//" "
1261 end do
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))
1268 end do
1269
1270 call io_close(iunit)
1271
1273 end subroutine ions_write_poscar
1274
1275 ! ---------------------------------------------------------
1276 subroutine ions_read_xyz(this, fname, comment)
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
1281
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)
1287
1288 push_sub(ions_read_xyz)
1289
1290 iunit = io_open(trim(fname)//'.xyz', this%namespace, action='read', position='rewind')
1291
1292 ! --- First two lines: natoms + optional comment ------------------------
1293 read(iunit, '(i4)') this%natoms
1294 if (present(comment)) then
1295 read(iunit, *)
1296 else
1297 read(iunit, *)
1298 end if
1299
1300 ! --- Read atom lines ---------------------------------------------------
1301 do iatom = 1, this%natoms
1302
1303 ! Read entire line (safe even if long)
1304 read(iunit,'(A)', iostat=ios) line
1305 if (ios /= 0) then
1306 call io_close(iunit)
1307 message(1) = "Error reading XYZ atom line"
1308 call messages_fatal(1, namespace=this%namespace)
1309 end if
1310
1311 ! Extract fixed-length label (preserves spaces inside)
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)
1316 end if
1317 label = line(1:label_len)
1318
1319 ! Extract remaining part for coordinates
1320 if (len(line) > label_len) then
1321 coordstr = adjustl(line(label_len+1:))
1322 else
1323 call io_close(iunit)
1324 message(1) = "XYZ file: no coordinate data after label"
1325 call messages_fatal(1, namespace=this%namespace)
1326 end if
1327
1328 ! Parse coordinates in a flexible (list-directed) way
1329 read(coordstr, *, iostat=ios) (tmp(idir), idir=1, this%space%dim)
1330 if (ios /= 0) then
1331 call io_close(iunit)
1332 message(1) = "XYZ file: malformed coordinate fields"
1333 call messages_fatal(1, namespace=this%namespace)
1334 end if
1335
1336 ! Convert and store
1337 this%pos(:, iatom) = units_to_atomic(units_out%length_xyz_file, tmp)
1338
1339 end do
1340
1341 call io_close(iunit)
1342
1343 pop_sub(ions_read_xyz)
1344 end subroutine ions_read_xyz
1345
1346 ! ----------------------------------------------------------------
1349 subroutine ions_write_crystal(this, dir)
1350 class(ions_t), intent(in) :: this
1351 character(len=*), intent(in) :: dir
1352
1353 type(lattice_iterator_t) :: latt_iter
1354 real(real64) :: radius, pos(this%space%dim)
1355 integer :: iatom, icopy, iunit
1356
1357 push_sub(ions_write_crystal)
1358
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)
1361
1362 if (mpi_world%is_root()) then
1363
1364 iunit = io_open(trim(dir)//'/crystal.xyz', this%namespace, action='write')
1365
1366 write(iunit, '(i9)') this%natoms*latt_iter%n_cells
1367 write(iunit, '(a)') '#generated by Octopus'
1368
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
1373 end do
1374 end do
1375
1376 call io_close(iunit)
1377 end if
1378
1379 pop_sub(ions_write_crystal)
1380 end subroutine ions_write_crystal
1381
1382 ! ---------------------------------------------------------
1383 subroutine ions_write_bild_forces_file(this, dir, fname)
1384 class(ions_t), intent(in) :: this
1385 character(len=*), intent(in) :: dir, fname
1386
1387 integer :: iunit, iatom, idir
1388 real(real64) :: force(this%space%dim), center(this%space%dim)
1389 character(len=20) frmt
1390
1391 if (.not. mpi_world%is_root()) return
1392
1394
1395 call io_mkdir(dir, this%namespace)
1396 iunit = io_open(trim(dir)//'/'//trim(fname)//'.bild', this%namespace, action='write', &
1397 position='asis')
1398
1399 write(frmt,'(a,i0,a)')'(a,2(', this%space%dim,'f16.6,1x))'
1400
1401 write(iunit, '(a)')'.comment : force vectors in ['//trim(units_abbrev(units_out%force))//']'
1402 write(iunit, *)
1403 write(iunit, '(a)')'.color red'
1404 write(iunit, *)
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)
1412 write(iunit,*)
1413 end do
1414
1415 call io_close(iunit)
1416
1418 end subroutine ions_write_bild_forces_file
1419
1420 ! -----------------------------------------------------
1421 subroutine ions_write_vtk_geometry(this, filename, ascii)
1422 class(ions_t), intent(in) :: this
1423 character(len=*), intent(in) :: filename
1424 logical, optional, intent(in) :: ascii
1425
1426 integer :: iunit, iatom, ierr
1427 logical :: ascii_
1428 real(real64), allocatable :: data(:, :)
1429 integer, allocatable :: idata(:, :)
1430 character(len=MAX_PATH_LEN) :: fullname
1431
1432 push_sub(ions_write_vtk_geometry)
1433
1434 assert(this%space%dim == 3)
1435
1436 ascii_ = optional_default(ascii, .true.)
1437
1438 fullname = trim(filename)//".vtk"
1439
1440 iunit = io_open(trim(fullname), this%namespace, action='write')
1441
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)
1445
1446 if (ascii_) then
1447 write(iunit, '(1a)') 'ASCII'
1448 else
1449 write(iunit, '(1a)') 'BINARY'
1450 end if
1451
1452 write(iunit, '(1a)') 'DATASET POLYDATA'
1453
1454 write(iunit, '(a,i9,a)') 'POINTS ', this%natoms, ' double'
1455
1456 if (ascii_) then
1457 do iatom = 1, this%natoms
1458 write(iunit, '(3f12.6)') this%pos(1:3, iatom)
1459 end do
1460 else
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)
1465 end do
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)') ''
1471 end if
1472
1473 write(iunit, '(a,2i9)') 'VERTICES ', this%natoms, 2*this%natoms
1474
1475 if (ascii_) then
1476 do iatom = 1, this%natoms
1477 write(iunit, '(2i9)') 1, iatom - 1
1478 end do
1479 else
1480 call io_close(iunit)
1481 safe_allocate(idata(1:2, 1:this%natoms))
1482 do iatom = 1, this%natoms
1483 idata(1, iatom) = 1
1484 idata(2, iatom) = iatom - 1
1485 end do
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)') ''
1491 end if
1492
1493 write(iunit, '(a,i9)') 'POINT_DATA', this%natoms
1494 write(iunit, '(a)') 'SCALARS element integer'
1495 write(iunit, '(a)') 'LOOKUP_TABLE default'
1496
1497 if (ascii_) then
1498 do iatom = 1, this%natoms
1499 write(iunit, '(i9)') nint(this%atom(iatom)%species%get_z())
1500 end do
1501 else
1502 call io_close(iunit)
1503
1504 safe_allocate(idata(1:this%natoms, 1))
1505
1506 do iatom = 1, this%natoms
1507 idata(iatom, 1) = nint(this%atom(iatom)%species%get_z())
1508 end do
1509
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())
1512
1513 safe_deallocate_a(idata)
1514
1515 iunit = io_open(trim(fullname), this%namespace, action='write', position = 'append')
1516 write(iunit, '(1a)') ''
1517 end if
1518
1519 call io_close(iunit)
1520
1522 end subroutine ions_write_vtk_geometry
1523
1524 ! ---------------------------------------------------------
1525 function ions_current(this) result(current)
1526 class(ions_t), intent(in) :: this
1527 real(real64) :: current(this%space%dim)
1528
1529 integer :: iatom
1530
1531 push_sub(ions_current)
1532
1533 current = m_zero
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)
1536 end do
1537
1538 if (this%atoms_dist%parallel) then
1539 call comm_allreduce(this%atoms_dist%mpi_grp, current, dim=this%space%dim)
1540 end if
1541
1542 pop_sub(ions_current)
1543 end function ions_current
1544
1545 ! ---------------------------------------------------------
1546 function ions_abs_current(this) result(abs_current)
1547 class(ions_t), intent(in) :: this
1548 real(real64) :: abs_current(this%space%dim)
1549
1550 integer :: iatom
1551
1552 push_sub(ions_abs_current)
1553
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))
1557 end do
1558
1559 if (this%atoms_dist%parallel) then
1560 call comm_allreduce(this%atoms_dist%mpi_grp, abs_current, dim=this%space%dim)
1561 end if
1562
1563 pop_sub(ions_abs_current)
1564 end function ions_abs_current
1565
1566 ! ---------------------------------------------------------
1567 subroutine ions_finalize(ions)
1568 type(ions_t), intent(inout) :: ions
1569
1570 integer :: i
1571
1572 push_sub(ions_finalize)
1573
1574 call distributed_end(ions%atoms_dist)
1575
1576 call ion_interaction_end(ions%ion_interaction)
1577
1578 safe_deallocate_a(ions%atom)
1579 ions%natoms=0
1580
1581 if(allocated(ions%species)) then
1582 do i = 1, ions%nspecies
1583 !SAFE_ DEALLOCATE_P(ions%species(i)%s)
1584 if(associated(ions%species(i)%s)) deallocate(ions%species(i)%s)
1585 end do
1586 deallocate(ions%species)
1587 end if
1588 ions%nspecies=0
1589
1590 safe_deallocate_a(ions%pos_displacements)
1591 safe_deallocate_a(ions%vel_displacements)
1592
1593 call charged_particles_end(ions)
1594
1595 safe_deallocate_a(ions%map_symm_atoms)
1596 safe_deallocate_a(ions%inv_map_symm_atoms)
1597
1598
1599 pop_sub(ions_finalize)
1600 end subroutine ions_finalize
1601
1602 !-------------------------------------------------------------------
1605 subroutine ions_update_lattice_vectors(ions, latt, symmetrize)
1606 class(ions_t), intent(inout) :: ions
1607 type(lattice_vectors_t), intent(in) :: latt
1608 logical, intent(in) :: symmetrize
1609
1611
1612 ! Update the lattice vectors
1613 call ions%latt%update(latt%rlattice)
1614
1615 ! Regenerate symmetries in Cartesian space
1616 if (symmetrize) then
1617 call symmetries_update_lattice_vectors(ions%symm, latt, ions%space%dim)
1618 end if
1619
1621 end subroutine ions_update_lattice_vectors
1622
1623 !-------------------------------------------------------------------
1630 class(ions_t), intent(inout) :: ions
1631
1632 integer :: iatom, iop, iatom_symm, dim4symms
1633 real(real64) :: ratom(ions%space%dim)
1634
1636
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))
1639
1640 ! In the 4D case, symmetries are only defined in 3D, see symmetries.F90
1641 dim4symms = min(3, ions%space%dim)
1642 ratom = m_zero
1643
1644 do iop = 1, ions%symm%nops
1645 do iatom = 1, ions%natoms
1646 !We find the atom that correspond to this one, once symmetry is applied
1647 ratom(1:dim4symms) = symm_op_apply_inv_cart(ions%symm%ops(iop), ions%pos(:, iatom))
1648
1649 ratom(:) = ions%latt%fold_into_cell(ratom(:))
1650
1651 ! find iatom_symm
1652 do iatom_symm = 1, ions%natoms
1653 if (all(abs(ratom(:) - ions%pos(:, iatom_symm)) < symprec)) exit
1654 end do
1655
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)
1660 end if
1661
1662 ions%map_symm_atoms(iatom, iop) = iatom_symm
1663 ions%inv_map_symm_atoms(iatom_symm, iop) = iatom
1664 end do
1665 end do
1666
1667 ! Also build a map for non-symmorphic operations
1668 do iop = 1, ions%symm%nops_nonsymmorphic
1669 do iatom = 1, ions%natoms
1670 !We find the atom that correspond to this one, once symmetry is applied
1671 ratom(1:dim4symms) = symm_op_apply_inv_cart(ions%symm%non_symmorphic_ops(iop), ions%pos(:, iatom))
1672
1673 ratom(:) = ions%latt%fold_into_cell(ratom(:))
1674
1675 ! find iatom_symm
1676 do iatom_symm = 1, ions%natoms
1677 if (all(abs(ratom(:) - ions%pos(:, iatom_symm)) < symprec)) exit
1678 end do
1679
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)
1684 end if
1685
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
1688 end do
1689 end do
1690
1691
1694
1695 !-------------------------------------------------------------------
1699 subroutine ions_symmetrize_atomic_coord(ions)
1700 class(ions_t), intent(inout) :: ions
1701
1702 integer :: iatom, iop, iatom_sym
1703 real(real64) :: ratom(ions%space%dim)
1704 real(real64), allocatable :: new_pos(:,:)
1705
1707
1708 safe_allocate(new_pos(1:ions%space%dim, 1:ions%natoms))
1709
1710 do iatom = 1, ions%natoms
1711 new_pos(:, iatom) = m_zero
1712
1713 ! Symmorphic operations
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
1719 end do
1720
1721 ! Non-symmorphic operations
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
1727 end do
1728
1729 new_pos(:, iatom) = new_pos(:, iatom) / (ions%symm%nops + ions%symm%nops_nonsymmorphic)
1730 end do
1731
1732
1734 end subroutine ions_symmetrize_atomic_coord
1735
1737 subroutine ions_print_spacegroup(ions)
1738 class(ions_t), intent(inout) :: ions
1739
1740 type(spglibdataset) :: spg_dataset
1741 character(len=11) :: symbol
1742 integer, allocatable :: site_type(:)
1743 integer :: space_group, ia
1744
1745 if(.not. ions%space%is_periodic()) return
1746
1747 push_sub(ions_print_spacegroup)
1748
1749 safe_allocate(site_type(1:ions%natoms))
1750 do ia = 1, ions%natoms
1751 site_type(ia) = ions%atom(ia)%species%get_index()
1752 end do
1753
1754 spg_dataset = symmetries_get_spg_dataset(ions%namespace, ions%space, ions%latt, ions%natoms, ions%pos, site_type)
1755
1756 safe_deallocate_a(site_type)
1757
1758 if (spg_dataset%spglib_error /= 0) then
1759 pop_sub(ions_print_spacegroup)
1760 return
1761 end if
1762
1763 space_group = spg_dataset%spacegroup_number
1764 symbol = spg_dataset%international_symbol
1765
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)
1769
1770 pop_sub(ions_print_spacegroup)
1771 end subroutine ions_print_spacegroup
1772
1773end module ions_oct_m
1774
1775!! Local Variables:
1776!! mode: f90
1777!! coding: utf-8
1778!! End:
--------------— axpy ---------------— Constant times a vector plus a vector.
Definition: blas.F90:178
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)
Definition: ions.F90:1214
subroutine, public atom_init(this, dim, label, species)
Definition: atom.F90:151
subroutine, public atom_get_species(this, species)
Definition: atom.F90:261
subroutine, public atom_set_species(this, species)
Definition: atom.F90:249
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:120
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
Definition: global.F90:209
real(real64), parameter, public m_zero
Definition: global.F90:191
character(len= *), parameter, public static_dir
Definition: global.F90:266
real(real64), parameter, public p_kb
Boltzmann constant in Ha/K.
Definition: global.F90:228
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
Definition: io.F90:116
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:361
subroutine, public ion_interaction_init_parallelization(this, natoms, mc)
subroutine, public ion_interaction_init(this, namespace, space, natoms)
subroutine ions_init_interaction(this, interaction)
Definition: ions.F90:681
subroutine ions_update_lattice_vectors(ions, latt, symmetrize)
Regenerate the ions information after update of the lattice vectors.
Definition: ions.F90:1701
class(ions_t) function, pointer ions_constructor(namespace, print_info, latt_inp)
Definition: ions.F90:241
subroutine ions_fold_atoms_into_cell(this)
Definition: ions.F90:999
subroutine ions_single_mode_displacements(ions, mode, amplitude_pos, amplitude_vel)
apply initial displacements and velocities corresponding to a given phonon mode
Definition: ions.F90:857
real(real64) function, dimension(this%space%dim) ions_global_force(this, time)
Definition: ions.F90:1256
subroutine ions_copy(ions_out, ions_in)
Definition: ions.F90:635
real(real64) function, dimension(this%space%dim) ions_current(this)
Definition: ions.F90:1621
subroutine ions_update_quantity(this, label)
Definition: ions.F90:952
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,...
Definition: ions.F90:1725
subroutine ions_write_xyz(this, fname, append, comment, reduce_coordinates)
Definition: ions.F90:1273
subroutine ions_init_species(ions, factory, print_info)
Definition: ions.F90:523
subroutine ions_finalize(ions)
Definition: ions.F90:1663
subroutine ions_initialize(this)
Definition: ions.F90:938
real(real64) function, dimension(this%space%dim) ions_abs_current(this)
Definition: ions.F90:1642
subroutine ions_read_xyz(this, fname, comment)
Definition: ions.F90:1372
subroutine ions_write_vtk_geometry(this, filename, ascii)
Definition: ions.F90:1517
subroutine ions_rotate(this, from, from2, to)
Definition: ions.F90:1142
subroutine ions_translate(this, xx)
Definition: ions.F90:1126
subroutine ions_symmetrize_atomic_coord(ions)
Symmetrizes atomic coordinates by applying all symmetries.
Definition: ions.F90:1795
subroutine ions_print_spacegroup(ions)
Prints the spacegroup of the system for periodic systems.
Definition: ions.F90:1833
real(real64) function ions_min_distance(this, real_atoms_only)
Definition: ions.F90:1014
subroutine ions_init_random_displacements(ions, T)
create random displacements for positions and velocities
Definition: ions.F90:704
real(real64) function, dimension(this%space%dim) ions_dipole(this, mask)
Definition: ions.F90:1104
subroutine ions_write_bild_forces_file(this, dir, fname)
Definition: ions.F90:1479
subroutine ions_write_crystal(this, dir)
This subroutine creates a crystal by replicating the geometry and writes the result to dir.
Definition: ions.F90:1445
logical function ions_has_time_dependent_species(this)
Definition: ions.F90:1077
subroutine ions_partition(this, mc)
Definition: ions.F90:667
subroutine ions_init_interaction_as_partner(partner, interaction)
Definition: ions.F90:968
subroutine ions_copy_quantities_to_interaction(partner, interaction)
Definition: ions.F90:983
subroutine ions_write_poscar(this, fname, comment)
Writes the positions of the ions in POSCAR format.
Definition: ions.F90:1326
real(real64) function ions_val_charge(this, mask)
Definition: ions.F90:1088
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
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.
Definition: tdfunction.F90:220
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
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.
Definition: multicomm.F90:208
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...
Definition: species.F90:147
Class describing a Wigner distribution for sampling initial conditions in multi-trajectory Ehrenfest ...
int true(void)