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 call mpi_grp_init(this%grp, mc%master_comm)
582
583 pop_sub(ions_partition)
584 end subroutine ions_partition
585
586 ! ---------------------------------------------------------
587 subroutine ions_init_interaction(this, interaction)
588 class(ions_t), target, intent(inout) :: this
589 class(interaction_t), intent(inout) :: interaction
590
591 push_sub(ions_init_interaction)
592
593 select type (interaction)
594 class default
595 call charged_particles_init_interaction(this, interaction)
596 end select
597
598 pop_sub(ions_init_interaction)
599 end subroutine ions_init_interaction
600
601 ! ---------------------------------------------------------
609 !
610 subroutine ions_init_random_displacements(ions, T)
611 class(ions_t), intent(inout) :: ions
612 real(real64), intent(in) :: T
613
614
615 integer(int64) :: seed
616
617 type(phonon_modes_t) :: phonons
618 type(wigner_distribution_t) :: wigner
619
620 real(real64), allocatable :: sigmaP(:)
621 real(real64), allocatable :: muP(:)
622 real(real64), allocatable :: sigmaQ(:)
623 real(real64), allocatable :: muQ(:)
624 real(real64), allocatable :: genP(:)
625 real(real64), allocatable :: genQ(:)
626
627
628 real(real64) :: beta_half
629 real(real64), allocatable :: pos_displacements(:), vel_displacements(:)
630 real(real64), allocatable :: l_m(:), prefactor(:)
631
632 integer :: imode, idim, iatom, i, iline
633 integer :: num_real_modes
635 real(real64), parameter :: low_temperature_tolerance = 1.0e-6_real64
636
637
639
640 ! create a unique seed for each replica, based on a hash of the full (unique) namespace string.
641 seed = ions%namespace%get_hash32()
642
643 message(1) = "Create initial random displacements for ions."
644 call messages_info(1, namespace = ions%namespace)
645
646 write(message(1), '("namespace = ",A)') ions%namespace%get()
647 write(message(2), '("seed = ",I0)') seed
648 call messages_info(2, namespace = ions%namespace, debug_only=.true.)
649
650
651 ! initialize the phonons (read info from file)
652 call phonons%init(ions%namespace, ions%space%dim, ions%natoms, ions%space%periodic_dim > 0)
653
654 num_real_modes = phonons%num_modes
655
656 if (num_real_modes == 0) then
657
658 ions%pos_displacements = m_zero
659 ions%vel_displacements = m_zero
660
662 return
663 end if
664
665
666 ! initialize wigner distribution for phonons%num_modes modes and with given seed
667 call wigner%init(num_real_modes, seed)
668
669 ! set up widths and shifts for the Wigner function
670 safe_allocate(sigmap(1:num_real_modes))
671 safe_allocate(sigmaq(1:num_real_modes))
672 safe_allocate(mup(1:num_real_modes))
673 safe_allocate(muq(1:num_real_modes))
674 safe_allocate(genp(1:num_real_modes))
675 safe_allocate(genq(1:num_real_modes))
676
677 safe_allocate(pos_displacements(1:(ions%space%dim*ions%natoms)))
678 safe_allocate(vel_displacements(1:(ions%space%dim*ions%natoms)))
679 safe_allocate(prefactor(1:(ions%space%dim*ions%natoms)))
680
681 safe_allocate(l_m(1:num_real_modes))
683
684 if (t < low_temperature_tolerance) then
685 ! TODO: Implement proper low T limit
686 sigmap = 1.0_real64/2.0_real64
687 sigmaq = 1.0_real64/2.0_real64
688 else
689 beta_half = m_one / (2 * p_kb * t)
690 sigmap = 1.0_real64/sqrt((2.0_real64 * tanh(beta_half * phonons%frequencies)))
691 sigmaq = 1.0_real64/sqrt((2.0_real64 * tanh(beta_half * phonons%frequencies)))
692 end if
693 mup = m_zero
694 muq = m_zero
695
696 l_m = sqrt(2.0_real64/(unit_amu%factor * phonons%frequencies)) * phonons%alpha
697
698 i = 1
699 do iatom = 1, ions%natoms
700 do idim = 1, ions%space%dim
701 prefactor(i) = 1.0_real64/sqrt(ions%mass(iatom)/unit_amu%factor * phonons%num_super)
702 i=i+1
703 end do
704 end do
706 ! get generalized mode coords from Wigner function
707 genq = wigner%get(sigmaq, muq, wigner_q) * l_m
708 genp = wigner%get(sigmap, mup, wigner_p) / (l_m * unit_amu%factor) ! TODO: Check!!
709
710 ! construct initial displacements and velocities for ions (velocities must be applied later)
711
712 ions%pos_displacements = m_zero
713 ions%vel_displacements = m_zero
714
715 ! implement eq.(9) from Kevins paper for periodic systems
716
717 do imode = 1, num_real_modes
718
719 pos_displacements = prefactor(:) * phonons%eigenvectors(:, imode) * genq(imode)
720 vel_displacements = prefactor(:) * phonons%eigenvectors(:, imode) * genp(imode)
721
722 iline = 1
723 write(message(iline), '("Displacements for mode ",I4)') imode
724 do iatom=1, ions%natoms
725 iline = iline+1
726 write(message(iline), '(3E15.5)') pos_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
727 end do
728 iline = iline+1
729 write(message(iline), '("Velocities for mode ",I4)') imode
730 do iatom=1, ions%natoms
731 iline = iline+1
732 write(message(iline), '(3E15.5)') vel_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
733 end do
734
735 call messages_info(iline, namespace=ions%namespace)
736
737 call blas_axpy(phonons%dim, 1.0_real64, pos_displacements(1), 1, ions%pos_displacements(1,1), 1)
738 call blas_axpy(phonons%dim, 1.0_real64, vel_displacements(1), 1, ions%vel_displacements(1,1), 1)
739
740 end do
741
742 safe_deallocate_a(sigmap)
743 safe_deallocate_a(sigmaq)
744 safe_deallocate_a(mup)
745 safe_deallocate_a(muq)
746 safe_deallocate_a(genp)
747 safe_deallocate_a(genq)
748 safe_deallocate_a(l_m)
749
750 safe_deallocate_a(prefactor)
751 safe_deallocate_a(pos_displacements)
752 safe_deallocate_a(vel_displacements)
753
755
756 end subroutine ions_init_random_displacements
757
758
762 !
763 subroutine ions_single_mode_displacements(ions, mode, amplitude_pos, amplitude_vel)
764 class(ions_t), intent(inout) :: ions
765 integer, intent(in) :: mode
766 real(real64), optional, intent(in) :: amplitude_pos
767 real(real64), optional, intent(in) :: amplitude_vel
768
769 type(phonon_modes_t) :: phonons
770 integer :: num_real_modes, i, iatom, idim
771 real(real64) :: l_m, genP, genQ
772
773 real(real64), allocatable :: pos_displacements(:), vel_displacements(:), prefactor(:)
774
775
777
778 write(message(1), '("Create displacements for mode ", I4)') mode
779 call messages_info(1, namespace=ions%namespace)
780
781 ! initialize the phonons (read info from file)
782 call phonons%init(ions%namespace, ions%space%dim, ions%natoms, ions%space%periodic_dim > 0)
783
784 num_real_modes = phonons%num_modes
785
786 if (mode > num_real_modes) then
787 write(message(1), '("Requested mode ",I0," exceeds number of available modes (",I0,")")') mode, num_real_modes
788 call messages_fatal(1, namespace=ions%namespace)
789
791 return
792 end if
793
794 ions%pos_displacements = m_zero
795 ions%vel_displacements = m_zero
796
797 l_m = sqrt(2.0_real64/(unit_amu%factor * phonons%frequencies(mode))) * phonons%alpha(mode)
798
799 safe_allocate(pos_displacements(1:(ions%space%dim*ions%natoms)))
800 safe_allocate(vel_displacements(1:(ions%space%dim*ions%natoms)))
801 safe_allocate(prefactor(1:(ions%space%dim*ions%natoms)))
802
803 i = 1
804 do iatom = 1, ions%natoms
805 do idim=1, ions%space%dim
806 prefactor(i) = 1.0_real64/sqrt(ions%mass(iatom)/unit_amu%factor * phonons%num_super)
807 i=i+1
808 end do
809 end do
810
811 ! get generalized mode coords from Wigner function
812 genq = amplitude_pos * l_m
813 genp = amplitude_vel / (l_m * unit_amu%factor)! TODO: Check!!
814
815 pos_displacements = prefactor(:) * phonons%eigenvectors(:, mode) * genq
816 vel_displacements = prefactor(:) * phonons%eigenvectors(:, mode) * genp
817
818 write(message(1), '("Displacements for mode ",I4)') mode
819 call messages_info(1, namespace=ions%namespace)
820 do iatom=1, ions%natoms
821 write(message(1), '(3E15.5)') pos_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
822 call messages_info(1, namespace=ions%namespace)
823 end do
824 write(message(1), '("Velocities for mode ",I4)') mode
825 call messages_info(1, namespace=ions%namespace)
826 do iatom=1, ions%natoms
827 write(message(1), '(3E15.5)') vel_displacements((iatom-1)*ions%space%dim+1:(iatom-1)*ions%space%dim+ions%space%dim)
828 call messages_info(1, namespace=ions%namespace)
829 end do
830
831
832 call blas_axpy(phonons%dim, 1.0_real64, pos_displacements(1), 1, ions%pos_displacements(1,1), 1)
833 call blas_axpy(phonons%dim, 1.0_real64, vel_displacements(1), 1, ions%vel_displacements(1,1), 1)
834
835 safe_deallocate_a(pos_displacements)
836 safe_deallocate_a(vel_displacements)
837 safe_deallocate_a(prefactor)
838
840
841 end subroutine ions_single_mode_displacements
842
843 ! ---------------------------------------------------------
844 subroutine ions_initialize(this)
845 class(ions_t), intent(inout) :: this
846
847 push_sub(ions_initialize)
848
849 ! At this point, we need to apply initial velocities from the random displacements.
850 if(allocated(this%vel_displacements)) then
851 this%vel = this%vel + this%vel_displacements
852 end if
853
854 pop_sub(ions_initialize)
855 end subroutine ions_initialize
856
857 ! ---------------------------------------------------------
858 subroutine ions_update_quantity(this, label)
859 class(ions_t), intent(inout) :: this
860 character(len=*), intent(in) :: label
861
862 push_sub(ions_update_quantity)
863
864 select case (label)
865 case default
866 ! Other quantities should be handled by the parent class
867 call charged_particles_update_quantity(this, label)
868 end select
869
870 pop_sub(ions_update_quantity)
871 end subroutine ions_update_quantity
872
873 ! ---------------------------------------------------------
874 subroutine ions_init_interaction_as_partner(partner, interaction)
875 class(ions_t), intent(in) :: partner
876 class(interaction_surrogate_t), intent(inout) :: interaction
877
879
880 select type (interaction)
881 class default
882 call charged_particles_init_interaction_as_partner(partner, interaction)
883 end select
884
887
888 ! ---------------------------------------------------------
889 subroutine ions_copy_quantities_to_interaction(partner, interaction)
890 class(ions_t), intent(inout) :: partner
891 class(interaction_surrogate_t), intent(inout) :: interaction
892
894
895 select type (interaction)
896 class default
897 ! Other interactions should be handled by the parent class
898 call charged_particles_copy_quantities_to_interaction(partner, interaction)
899 end select
900
903
904 ! ---------------------------------------------------------
905 subroutine ions_fold_atoms_into_cell(this)
906 class(ions_t), intent(inout) :: this
907
908 integer :: iatom
909
911
912 do iatom = 1, this%natoms
913 this%pos(:, iatom) = this%latt%fold_into_cell(this%pos(:, iatom))
914 end do
915
917 end subroutine ions_fold_atoms_into_cell
918
919 ! ---------------------------------------------------------
920 real(real64) function ions_min_distance(this, real_atoms_only) result(rmin)
921 class(ions_t), intent(in) :: this
922 logical, optional, intent(in) :: real_atoms_only
923
924 integer :: iatom, jatom, idir
925 real(real64) :: xx(this%space%dim)
926 logical :: real_atoms_only_
927 class(species_t), pointer :: species
928
929 if (this%natoms == 1 .and. .not. this%space%is_periodic()) then
930 rmin = m_huge
931 return
932 end if
933
934 push_sub(ions_min_distance)
935
936 real_atoms_only_ = optional_default(real_atoms_only, .false.)
937
938 ! Without this line, valgrind complains about a conditional jump on uninitialized variable,
939 ! as atom_get_species has an intent(out) that causes a call to the finalizer (with Ifort)
940 nullify(species)
941
942 rmin = huge(rmin)
943 do iatom = 1, this%natoms
944 call atom_get_species(this%atom(iatom), species)
945 select type(species)
946 class is(jellium_t)
947 if (real_atoms_only_) cycle
948 end select
949 do jatom = iatom + 1, this%natoms
950 call atom_get_species(this%atom(jatom), species)
951 select type(species)
952 class is(jellium_t)
953 if (real_atoms_only_) cycle
954 end select
955 xx = abs(this%pos(:, iatom) - this%pos(:, jatom))
956 if (this%space%is_periodic()) then
957 xx = this%latt%cart_to_red(xx)
958 do idir = 1, this%space%periodic_dim
959 xx(idir) = xx(idir) - floor(xx(idir) + m_half)
960 end do
961 xx = this%latt%red_to_cart(xx)
962 end if
963 rmin = min(norm2(xx), rmin)
964 end do
965 end do
966
967 if (.not. (this%only_user_def .and. real_atoms_only_)) then
968 ! what if the nearest neighbors are periodic images?
969 do idir = 1, this%space%periodic_dim
970 rmin = min(rmin, norm2(this%latt%rlattice(:,idir)))
971 end do
972 end if
973
974 ! To avoid numerical instabilities, we round this to 6 digits only
975 if(rmin < huge(rmin)/1e6_real64) then
976 rmin = anint(rmin*1e6_real64)*1.0e-6_real64
977 end if
978
979 pop_sub(ions_min_distance)
980 end function ions_min_distance
981
982 ! ---------------------------------------------------------
983 logical function ions_has_time_dependent_species(this) result(time_dependent)
984 class(ions_t), intent(in) :: this
985
987
988 time_dependent = this%species_time_dependent
989
992
993 ! ---------------------------------------------------------
994 real(real64) function ions_val_charge(this, mask) result(val_charge)
995 class(ions_t), intent(in) :: this
996 logical, optional, intent(in) :: mask(:)
997
998 push_sub(ions_val_charge)
999
1000 if (present(mask)) then
1001 val_charge = -sum(this%charge, mask=mask)
1002 else
1003 val_charge = -sum(this%charge)
1004 end if
1005
1006 pop_sub(ions_val_charge)
1007 end function ions_val_charge
1008
1009 ! ---------------------------------------------------------
1010 function ions_dipole(this, mask) result(dipole)
1011 class(ions_t), intent(in) :: this
1012 logical, optional, intent(in) :: mask(:)
1013 real(real64) :: dipole(this%space%dim)
1014
1015 integer :: ia
1016
1017 push_sub(ions_dipole)
1018
1019 dipole = m_zero
1020 do ia = 1, this%natoms
1021 if (present(mask)) then
1022 if (.not. mask(ia)) cycle
1023 end if
1024 dipole = dipole + this%charge(ia)*this%pos(:, ia)
1025 end do
1026 dipole = p_proton_charge*dipole
1027
1028 pop_sub(ions_dipole)
1029 end function ions_dipole
1030
1031 ! ---------------------------------------------------------
1032 subroutine ions_translate(this, xx)
1033 class(ions_t), intent(inout) :: this
1034 real(real64), intent(in) :: xx(this%space%dim)
1035
1036 integer :: iatom
1037
1038 push_sub(ions_translate)
1039
1040 do iatom = 1, this%natoms
1041 this%pos(:, iatom) = this%pos(:, iatom) - xx
1042 end do
1043
1044 pop_sub(ions_translate)
1045 end subroutine ions_translate
1046
1047 ! ---------------------------------------------------------
1048 subroutine ions_rotate(this, from, from2, to)
1049 class(ions_t), intent(inout) :: this
1050 real(real64), intent(in) :: from(this%space%dim)
1051 real(real64), intent(in) :: from2(this%space%dim)
1052 real(real64), intent(in) :: to(this%space%dim)
1053
1054 integer :: iatom
1055 real(real64) :: m1(3, 3), m2(3, 3)
1056 real(real64) :: m3(3, 3), f2(3), per(3)
1057 real(real64) :: alpha, r
1058
1059 push_sub(ions_rotate)
1060
1061 if (this%space%dim /= 3) then
1062 call messages_not_implemented("ions_rotate in other than 3 dimensions", namespace=this%namespace)
1063 end if
1064
1065 ! initialize matrices
1066 m1 = diagonal_matrix(3, m_one)
1067
1068 ! rotate the to-axis to the z-axis
1069 if (abs(to(2)) > 1d-150) then
1070 alpha = atan2(to(2), to(1))
1071 call rotate(m1, alpha, 3)
1072 end if
1073 alpha = atan2(norm2(to(1:2)), to(3))
1074 call rotate(m1, -alpha, 2)
1075
1076 ! get perpendicular to z and from
1077 f2 = matmul(m1, from)
1078 per(1) = -f2(2)
1079 per(2) = f2(1)
1080 per(3) = m_zero
1081 r = norm2(per)
1082 if (r > m_zero) then
1083 per = per/r
1084 else
1085 per(2) = m_one
1086 end if
1087
1088 ! rotate perpendicular axis to the y-axis
1089 m2 = diagonal_matrix(3, m_one)
1090 alpha = atan2(per(1), per(2))
1091 call rotate(m2, -alpha, 3)
1092
1093 ! rotate from => to (around the y-axis)
1094 m3 = diagonal_matrix(3, m_one)
1095 alpha = acos(min(sum(from*to), m_one))
1096 call rotate(m3, -alpha, 2)
1097
1098 ! join matrices
1099 m2 = matmul(transpose(m2), matmul(m3, m2))
1100
1101 ! rotate around the z-axis to get the second axis
1102 per = matmul(m2, matmul(m1, from2))
1103 alpha = atan2(per(1), per(2))
1104 call rotate(m2, -alpha, 3) ! second axis is now y
1106 ! get combined transformation
1107 m1 = matmul(transpose(m1), matmul(m2, m1))
1108
1109 ! now transform the coordinates
1110 ! it is written in this way to avoid what I consider a bug in the Intel compiler
1111 do iatom = 1, this%natoms
1112 f2 = this%pos(:, iatom)
1113 this%pos(:, iatom) = matmul(m1, f2)
1114 end do
1115
1116 pop_sub(ions_rotate)
1117 contains
1118
1119 ! ---------------------------------------------------------
1120 subroutine rotate(m, angle, dir)
1121 real(real64), intent(inout) :: m(3, 3)
1122 real(real64), intent(in) :: angle
1123 integer, intent(in) :: dir
1124
1125 real(real64) :: aux(3, 3), ca, sa
1126
1128
1129 ca = cos(angle)
1130 sa = sin(angle)
1131
1132 aux = m_zero
1133 select case (dir)
1134 case (1)
1135 aux(1, 1) = m_one
1136 aux(2, 2) = ca
1137 aux(3, 3) = ca
1138 aux(2, 3) = sa
1139 aux(3, 2) = -sa
1140 case (2)
1141 aux(2, 2) = m_one
1142 aux(1, 1) = ca
1143 aux(3, 3) = ca
1144 aux(1, 3) = sa
1145 aux(3, 1) = -sa
1146 case (3)
1147 aux(3, 3) = m_one
1148 aux(1, 1) = ca
1149 aux(2, 2) = ca
1150 aux(1, 2) = sa
1151 aux(2, 1) = -sa
1152 end select
1153
1154 m = matmul(aux, m)
1155
1156 pop_sub(ions_rotate.rotate)
1157 end subroutine rotate
1158
1159 end subroutine ions_rotate
1160
1161 ! ---------------------------------------------------------
1162 function ions_global_force(this, time) result(force)
1163 class(ions_t), intent(in) :: this
1164 real(real64), intent(in) :: time
1165 real(real64) :: force(this%space%dim)
1166
1167 push_sub(ions_global_force)
1168
1169 force = m_zero
1170
1171 if (this%apply_global_force) then
1172 force(1) = tdf(this%global_force_function, time)
1173 end if
1174
1175 pop_sub(ions_global_force)
1176 end function ions_global_force
1177
1178 ! ---------------------------------------------------------
1179 subroutine ions_write_xyz(this, fname, append, comment, reduce_coordinates)
1180 class(ions_t), intent(in) :: this
1181 character(len=*), intent(in) :: fname
1182 logical, optional, intent(in) :: append
1183 character(len=*), optional, intent(in) :: comment
1184 logical, optional, intent(in) :: reduce_coordinates
1185
1186 integer :: iatom, idim, iunit
1187 character(len=6) position
1188 character(len=19) :: frmt
1189 real(real64) :: red(this%space%dim)
1190
1191 if (.not. this%grp%is_root()) return
1192
1193 push_sub(ions_write_xyz)
1194
1195 position = 'asis'
1196 if (present(append)) then
1197 if (append) position = 'append'
1198 end if
1199 if(.not.optional_default(reduce_coordinates, .false.)) then
1200 iunit = io_open(trim(fname)//'.xyz', this%namespace, action='write', position=position)
1201 else
1202 iunit = io_open(trim(fname)//'.xyz_red', this%namespace, action='write', position=position)
1203 end if
1204
1205 write(iunit, '(i4)') this%natoms
1206 if (present(comment)) then
1207 write(iunit, '(1x,a)') comment
1208 else
1209 write(iunit, '(1x,a,a)') 'units: ', trim(units_abbrev(units_out%length_xyz_file))
1210 end if
1211
1212 write(unit=frmt, fmt="(a5,i2.2,a4,i2.2,a6)") "(6x,a", label_len, ",2x,", this%space%dim,"f15.9)"
1213 do iatom = 1, this%natoms
1214 if(.not.optional_default(reduce_coordinates, .false.)) then
1215 write(unit=iunit, fmt=frmt) this%atom(iatom)%label, &
1216 (units_from_atomic(units_out%length_xyz_file, this%pos(idim, iatom)), idim=1, this%space%dim)
1217 else
1218 red = this%latt%cart_to_red(this%pos(:, iatom))
1219 write(unit=iunit, fmt=frmt) this%atom(iatom)%label, (red(idim), idim=1, this%space%dim)
1220 end if
1221 end do
1222 call io_close(iunit)
1223
1224 pop_sub(ions_write_xyz)
1225 end subroutine ions_write_xyz
1226
1231 !
1232 subroutine ions_write_poscar(this, fname, comment)
1233 class(ions_t), intent(in) :: this
1234 character(len=*), optional, intent(in) :: fname
1235 character(len=*), optional, intent(in) :: comment
1236
1237 integer :: iatom, idim, iunit
1238 character(len=:), allocatable :: fname_, comment_
1239 character(len=10) :: format_string
1240
1241 push_sub(ions_write_poscar)
1242
1243 if (.not. this%grp%is_root()) then
1244 pop_sub(ions_write_poscar)
1245 return
1246 end if
1247
1248 comment_ = optional_default(comment, "")
1249 fname_ = optional_default(fname, "POSCAR")
1250
1251 iunit = io_open(trim(fname_), this%namespace, action='write')
1252
1253 write(iunit, '(A)') comment_ ! mandatory comment
1254 write(iunit, '("1.0")') ! scaling factor: we use 1 as lattice vectors are absolute
1255
1256 write(format_string, '("(",I1,A,")")') this%space%dim, 'F15.10'
1257 do idim=1, this%space%dim
1258 write(iunit, format_string) units_from_atomic(unit_angstrom, this%latt%rlattice(:,idim))
1259 end do
1260
1261 do iatom = 1, this%natoms
1262 write(iunit, '(A)', advance='NO') trim(this%atom(iatom)%label)//" "
1263 end do
1264 write(iunit, '("")')
1265 write(iunit, '(A)') repeat("1 ", this%natoms)
1266 write(iunit, '("direct")')
1267 write(format_string, '("(",I1,A,")")') this%space%dim, 'F15.10'
1268 do iatom=1, this%natoms
1269 write(iunit, format_string) this%latt%cart_to_red(this%pos(:, iatom))
1270 end do
1271
1272 call io_close(iunit)
1273
1275 end subroutine ions_write_poscar
1276
1277 ! ---------------------------------------------------------
1278 subroutine ions_read_xyz(this, fname, comment)
1279 use iso_fortran_env, only: real64
1280 class(ions_t), intent(inout) :: this
1281 character(len=*), intent(in) :: fname
1282 character(len=*), optional, intent(in) :: comment
1283
1284 integer :: iatom, idir, iunit, ios
1285 character(len=256) :: line
1286 character(len=LABEL_LEN) :: label
1287 character(len=:), allocatable :: coordstr
1288 real(real64) :: tmp(this%space%dim)
1289
1290 push_sub(ions_read_xyz)
1291
1292 iunit = io_open(trim(fname)//'.xyz', this%namespace, action='read', position='rewind')
1293
1294 ! --- First two lines: natoms + optional comment ------------------------
1295 read(iunit, '(i4)') this%natoms
1296 if (present(comment)) then
1297 read(iunit, *)
1298 else
1299 read(iunit, *)
1300 end if
1301
1302 ! --- Read atom lines ---------------------------------------------------
1303 do iatom = 1, this%natoms
1304
1305 ! Read entire line (safe even if long)
1306 read(iunit,'(A)', iostat=ios) line
1307 if (ios /= 0) then
1308 call io_close(iunit)
1309 message(1) = "Error reading XYZ atom line"
1310 call messages_fatal(1, namespace=this%namespace)
1311 end if
1312
1313 ! Extract fixed-length label (preserves spaces inside)
1314 if (len_trim(line) < label_len) then
1315 call io_close(iunit)
1316 message(1) = "XYZ file: atom line too short for label"
1317 call messages_fatal(1, namespace=this%namespace)
1318 end if
1319 label = line(1:label_len)
1320
1321 ! Extract remaining part for coordinates
1322 if (len(line) > label_len) then
1323 coordstr = adjustl(line(label_len+1:))
1324 else
1325 call io_close(iunit)
1326 message(1) = "XYZ file: no coordinate data after label"
1327 call messages_fatal(1, namespace=this%namespace)
1328 end if
1329
1330 ! Parse coordinates in a flexible (list-directed) way
1331 read(coordstr, *, iostat=ios) (tmp(idir), idir=1, this%space%dim)
1332 if (ios /= 0) then
1333 call io_close(iunit)
1334 message(1) = "XYZ file: malformed coordinate fields"
1335 call messages_fatal(1, namespace=this%namespace)
1336 end if
1337
1338 ! Convert and store
1339 this%pos(:, iatom) = units_to_atomic(units_out%length_xyz_file, tmp)
1340
1341 end do
1342
1343 call io_close(iunit)
1344
1345 pop_sub(ions_read_xyz)
1346 end subroutine ions_read_xyz
1347
1348 ! ----------------------------------------------------------------
1351 subroutine ions_write_crystal(this, dir)
1352 class(ions_t), intent(in) :: this
1353 character(len=*), intent(in) :: dir
1354
1355 type(lattice_iterator_t) :: latt_iter
1356 real(real64) :: radius, pos(this%space%dim)
1357 integer :: iatom, icopy, iunit
1358
1359 push_sub(ions_write_crystal)
1360
1361 radius = maxval(m_half*norm2(this%latt%rlattice(:,1:this%space%periodic_dim), dim=1))*(m_one + m_epsilon)
1362 latt_iter = lattice_iterator_t(this%latt, radius)
1363
1364 if (this%grp%is_root()) then
1365
1366 iunit = io_open(trim(dir)//'/crystal.xyz', this%namespace, action='write')
1367
1368 write(iunit, '(i9)') this%natoms*latt_iter%n_cells
1369 write(iunit, '(a)') '#generated by Octopus'
1370
1371 do iatom = 1, this%natoms
1372 do icopy = 1, latt_iter%n_cells
1373 pos = units_from_atomic(units_out%length, this%pos(:, iatom) + latt_iter%get(icopy))
1374 write(iunit, '(a, 99f12.6)') this%atom(iatom)%label, pos
1375 end do
1376 end do
1377
1378 call io_close(iunit)
1379 end if
1380
1381 pop_sub(ions_write_crystal)
1382 end subroutine ions_write_crystal
1383
1384 ! ---------------------------------------------------------
1385 subroutine ions_write_bild_forces_file(this, dir, fname)
1386 class(ions_t), intent(in) :: this
1387 character(len=*), intent(in) :: dir, fname
1388
1389 integer :: iunit, iatom, idir
1390 real(real64) :: force(this%space%dim), center(this%space%dim)
1391 character(len=20) frmt
1392
1393 if (.not. this%grp%is_root()) return
1394
1396
1397 call io_mkdir(dir, this%namespace)
1398 iunit = io_open(trim(dir)//'/'//trim(fname)//'.bild', this%namespace, action='write', &
1399 position='asis')
1400
1401 write(frmt,'(a,i0,a)')'(a,2(', this%space%dim,'f16.6,1x))'
1402
1403 write(iunit, '(a)')'.comment : force vectors in ['//trim(units_abbrev(units_out%force))//']'
1404 write(iunit, *)
1405 write(iunit, '(a)')'.color red'
1406 write(iunit, *)
1407 do iatom = 1, this%natoms
1408 center = units_from_atomic(units_out%length, this%pos(:, iatom))
1409 force = units_from_atomic(units_out%force, this%tot_force(:, iatom))
1410 write(iunit, '(a,1x,i4,1x,a2,1x,a6,1x,f10.6,a)') '.comment :', iatom, trim(this%atom(iatom)%label), &
1411 'force:', norm2(force), '['//trim(units_abbrev(units_out%force))//']'
1412 write(iunit,fmt=trim(frmt)) '.arrow', (center(idir), idir = 1, this%space%dim), &
1413 (center(idir) + force(idir), idir = 1, this%space%dim)
1414 write(iunit,*)
1415 end do
1416
1417 call io_close(iunit)
1418
1420 end subroutine ions_write_bild_forces_file
1421
1422 ! -----------------------------------------------------
1423 subroutine ions_write_vtk_geometry(this, filename, ascii)
1424 class(ions_t), intent(in) :: this
1425 character(len=*), intent(in) :: filename
1426 logical, optional, intent(in) :: ascii
1427
1428 integer :: iunit, iatom, ierr
1429 logical :: ascii_
1430 real(real64), allocatable :: data(:, :)
1431 integer, allocatable :: idata(:, :)
1432 character(len=MAX_PATH_LEN) :: fullname
1433
1434 push_sub(ions_write_vtk_geometry)
1435
1436 assert(this%space%dim == 3)
1437
1438 ascii_ = optional_default(ascii, .true.)
1439
1440 fullname = trim(filename)//".vtk"
1441
1442 iunit = io_open(trim(fullname), this%namespace, action='write')
1443
1444 write(iunit, '(1a)') '# vtk DataFile Version 3.0 '
1445 write(iunit, '(6a)') 'Generated by octopus ', trim(conf%version), ' - git: ', &
1446 trim(conf%git_commit), " configuration: ", trim(conf%config_time)
1447
1448 if (ascii_) then
1449 write(iunit, '(1a)') 'ASCII'
1450 else
1451 write(iunit, '(1a)') 'BINARY'
1452 end if
1453
1454 write(iunit, '(1a)') 'DATASET POLYDATA'
1455
1456 write(iunit, '(a,i9,a)') 'POINTS ', this%natoms, ' double'
1457
1458 if (ascii_) then
1459 do iatom = 1, this%natoms
1460 write(iunit, '(3f12.6)') this%pos(1:3, iatom)
1461 end do
1462 else
1463 call io_close(iunit)
1464 safe_allocate(data(1:3, 1:this%natoms))
1465 do iatom = 1, this%natoms
1466 data(1:3, iatom) = this%pos(1:3, iatom)
1467 end do
1468 call io_binary_write(io_workpath(fullname, this%namespace), i4_to_i8(3*this%natoms), data, &
1469 ierr, nohead = .true., fendian = io_binary_is_little_endian())
1470 safe_deallocate_a(data)
1471 iunit = io_open(trim(fullname), this%namespace, action='write', position = 'append')
1472 write(iunit, '(1a)') ''
1473 end if
1474
1475 write(iunit, '(a,2i9)') 'VERTICES ', this%natoms, 2*this%natoms
1476
1477 if (ascii_) then
1478 do iatom = 1, this%natoms
1479 write(iunit, '(2i9)') 1, iatom - 1
1480 end do
1481 else
1482 call io_close(iunit)
1483 safe_allocate(idata(1:2, 1:this%natoms))
1484 do iatom = 1, this%natoms
1485 idata(1, iatom) = 1
1486 idata(2, iatom) = iatom - 1
1487 end do
1488 call io_binary_write(io_workpath(fullname, this%namespace), i4_to_i8(2*this%natoms), idata, &
1489 ierr, nohead = .true., fendian = io_binary_is_little_endian())
1490 safe_deallocate_a(idata)
1491 iunit = io_open(trim(fullname), this%namespace, action='write', position = 'append')
1492 write(iunit, '(1a)') ''
1493 end if
1494
1495 write(iunit, '(a,i9)') 'POINT_DATA', this%natoms
1496 write(iunit, '(a)') 'SCALARS element integer'
1497 write(iunit, '(a)') 'LOOKUP_TABLE default'
1498
1499 if (ascii_) then
1500 do iatom = 1, this%natoms
1501 write(iunit, '(i9)') nint(this%atom(iatom)%species%get_z())
1502 end do
1503 else
1504 call io_close(iunit)
1505
1506 safe_allocate(idata(1:this%natoms, 1))
1507
1508 do iatom = 1, this%natoms
1509 idata(iatom, 1) = nint(this%atom(iatom)%species%get_z())
1510 end do
1511
1512 call io_binary_write(io_workpath(fullname, this%namespace), i4_to_i8(this%natoms), idata, &
1513 ierr, nohead = .true., fendian = io_binary_is_little_endian())
1514
1515 safe_deallocate_a(idata)
1516
1517 iunit = io_open(trim(fullname), this%namespace, action='write', position = 'append')
1518 write(iunit, '(1a)') ''
1519 end if
1520
1521 call io_close(iunit)
1522
1524 end subroutine ions_write_vtk_geometry
1525
1526 ! ---------------------------------------------------------
1527 function ions_current(this) result(current)
1528 class(ions_t), intent(in) :: this
1529 real(real64) :: current(this%space%dim)
1530
1531 integer :: iatom
1532
1533 push_sub(ions_current)
1534
1535 current = m_zero
1536 do iatom = this%atoms_dist%start, this%atoms_dist%end
1537 current = current + p_proton_charge*this%atom(iatom)%species%get_zval()*this%vel(:,iatom)
1538 end do
1539
1540 if (this%atoms_dist%parallel) then
1541 call comm_allreduce(this%atoms_dist%mpi_grp, current, dim=this%space%dim)
1542 end if
1543
1544 pop_sub(ions_current)
1545 end function ions_current
1546
1547 ! ---------------------------------------------------------
1548 function ions_abs_current(this) result(abs_current)
1549 class(ions_t), intent(in) :: this
1550 real(real64) :: abs_current(this%space%dim)
1551
1552 integer :: iatom
1553
1554 push_sub(ions_abs_current)
1555
1556 abs_current = m_zero
1557 do iatom = this%atoms_dist%start, this%atoms_dist%end
1558 abs_current = abs_current + abs(p_proton_charge*this%atom(iatom)%species%get_zval()*this%vel(:,iatom))
1559 end do
1560
1561 if (this%atoms_dist%parallel) then
1562 call comm_allreduce(this%atoms_dist%mpi_grp, abs_current, dim=this%space%dim)
1563 end if
1564
1565 pop_sub(ions_abs_current)
1566 end function ions_abs_current
1567
1568 ! ---------------------------------------------------------
1569 subroutine ions_finalize(ions)
1570 type(ions_t), intent(inout) :: ions
1571
1572 integer :: i
1573
1574 push_sub(ions_finalize)
1575
1576 call distributed_end(ions%atoms_dist)
1577
1578 call ion_interaction_end(ions%ion_interaction)
1579
1580 safe_deallocate_a(ions%atom)
1581 ions%natoms=0
1582
1583 if(allocated(ions%species)) then
1584 do i = 1, ions%nspecies
1585 !SAFE_ DEALLOCATE_P(ions%species(i)%s)
1586 if(associated(ions%species(i)%s)) deallocate(ions%species(i)%s)
1587 end do
1588 deallocate(ions%species)
1589 end if
1590 ions%nspecies=0
1591
1592 safe_deallocate_a(ions%pos_displacements)
1593 safe_deallocate_a(ions%vel_displacements)
1594
1595 call charged_particles_end(ions)
1596
1597 safe_deallocate_a(ions%map_symm_atoms)
1598 safe_deallocate_a(ions%inv_map_symm_atoms)
1599
1600
1601 pop_sub(ions_finalize)
1602 end subroutine ions_finalize
1603
1604 !-------------------------------------------------------------------
1607 subroutine ions_update_lattice_vectors(ions, latt, symmetrize)
1608 class(ions_t), intent(inout) :: ions
1609 type(lattice_vectors_t), intent(in) :: latt
1610 logical, intent(in) :: symmetrize
1611
1613
1614 ! Update the lattice vectors
1615 call ions%latt%update(latt%rlattice)
1616
1617 ! Regenerate symmetries in Cartesian space
1618 if (symmetrize) then
1619 call symmetries_update_lattice_vectors(ions%symm, latt, ions%space%dim)
1620 end if
1621
1623 end subroutine ions_update_lattice_vectors
1624
1625 !-------------------------------------------------------------------
1632 class(ions_t), intent(inout) :: ions
1633
1634 integer :: iatom, iop, iatom_symm, dim4symms
1635 real(real64) :: ratom(ions%space%dim)
1636
1638
1639 safe_allocate(ions%map_symm_atoms(1:ions%natoms, 1:ions%symm%nops + ions%symm%nops_nonsymmorphic))
1640 safe_allocate(ions%inv_map_symm_atoms(1:ions%natoms, 1:ions%symm%nops + ions%symm%nops_nonsymmorphic))
1641
1642 ! In the 4D case, symmetries are only defined in 3D, see symmetries.F90
1643 dim4symms = min(3, ions%space%dim)
1644 ratom = m_zero
1645
1646 do iop = 1, ions%symm%nops
1647 do iatom = 1, ions%natoms
1648 !We find the atom that correspond to this one, once symmetry is applied
1649 ratom(1:dim4symms) = symm_op_apply_inv_cart(ions%symm%ops(iop), ions%pos(:, iatom))
1650
1651 ratom(:) = ions%latt%fold_into_cell(ratom(:))
1652
1653 ! find iatom_symm
1654 do iatom_symm = 1, ions%natoms
1655 if (all(abs(ratom(:) - ions%pos(:, iatom_symm)) < symprec)) exit
1656 end do
1657
1658 if (iatom_symm > ions%natoms) then
1659 write(message(1),'(a,i6)') 'Internal error: could not find symmetric partner for atom number', iatom
1660 write(message(2),'(a,i3,a)') 'with symmetry operation number ', iop, '.'
1661 call messages_fatal(2, namespace=ions%namespace)
1662 end if
1663
1664 ions%map_symm_atoms(iatom, iop) = iatom_symm
1665 ions%inv_map_symm_atoms(iatom_symm, iop) = iatom
1666 end do
1667 end do
1668
1669 ! Also build a map for non-symmorphic operations
1670 do iop = 1, ions%symm%nops_nonsymmorphic
1671 do iatom = 1, ions%natoms
1672 !We find the atom that correspond to this one, once symmetry is applied
1673 ratom(1:dim4symms) = symm_op_apply_inv_cart(ions%symm%non_symmorphic_ops(iop), ions%pos(:, iatom))
1674
1675 ratom(:) = ions%latt%fold_into_cell(ratom(:))
1676
1677 ! find iatom_symm
1678 do iatom_symm = 1, ions%natoms
1679 if (all(abs(ratom(:) - ions%pos(:, iatom_symm)) < symprec)) exit
1680 end do
1681
1682 if (iatom_symm > ions%natoms) then
1683 write(message(1),'(a,i6)') 'Internal error: could not find symmetric partner for atom number', iatom
1684 write(message(2),'(a,i3,a)') 'with symmetry operation number ', iop, '.'
1685 call messages_fatal(2, namespace=ions%namespace)
1686 end if
1687
1688 ions%map_symm_atoms(iatom, ions%symm%nops+iop) = iatom_symm
1689 ions%inv_map_symm_atoms(iatom_symm, ions%symm%nops+iop) = iatom
1690 end do
1691 end do
1692
1693
1696
1697 !-------------------------------------------------------------------
1701 subroutine ions_symmetrize_atomic_coord(ions)
1702 class(ions_t), intent(inout) :: ions
1703
1704 integer :: iatom, iop, iatom_sym
1705 real(real64) :: ratom(ions%space%dim)
1706 real(real64), allocatable :: new_pos(:,:)
1707
1709
1710 safe_allocate(new_pos(1:ions%space%dim, 1:ions%natoms))
1711
1712 do iatom = 1, ions%natoms
1713 new_pos(:, iatom) = m_zero
1714
1715 ! Symmorphic operations
1716 do iop = 1, ions%symm%nops
1717 iatom_sym = ions%inv_map_symm_atoms(iatom, iop)
1718 ratom = symm_op_apply_inv_cart(ions%symm%ops(iop), ions%pos(:, iatom_sym))
1719 ratom = ions%latt%fold_into_cell(ratom)
1720 new_pos(:, iatom) = new_pos(:, iatom) + ratom
1721 end do
1722
1723 ! Non-symmorphic operations
1724 do iop = 1, ions%symm%nops_nonsymmorphic
1725 iatom_sym = ions%inv_map_symm_atoms(iatom, iop + ions%symm%nops)
1726 ratom = symm_op_apply_inv_cart(ions%symm%non_symmorphic_ops(iop), ions%pos(:, iatom_sym))
1727 ratom = ions%latt%fold_into_cell(ratom)
1728 new_pos(:, iatom) = new_pos(:, iatom) + ratom
1729 end do
1730
1731 new_pos(:, iatom) = new_pos(:, iatom) / (ions%symm%nops + ions%symm%nops_nonsymmorphic)
1732 end do
1733
1734
1736 end subroutine ions_symmetrize_atomic_coord
1737
1739 subroutine ions_print_spacegroup(ions)
1740 class(ions_t), intent(inout) :: ions
1741
1742 type(spglibdataset) :: spg_dataset
1743 character(len=11) :: symbol
1744 integer, allocatable :: site_type(:)
1745 integer :: space_group, ia
1746
1747 if(.not. ions%space%is_periodic()) return
1748
1749 push_sub(ions_print_spacegroup)
1750
1751 safe_allocate(site_type(1:ions%natoms))
1752 do ia = 1, ions%natoms
1753 site_type(ia) = ions%atom(ia)%species%get_index()
1754 end do
1755
1756 spg_dataset = symmetries_get_spg_dataset(ions%namespace, ions%space, ions%latt, ions%natoms, ions%pos, site_type)
1757
1758 safe_deallocate_a(site_type)
1759
1760 if (spg_dataset%spglib_error /= 0) then
1761 pop_sub(ions_print_spacegroup)
1762 return
1763 end if
1764
1765 space_group = spg_dataset%spacegroup_number
1766 symbol = spg_dataset%international_symbol
1767
1768 write(message(1),'(a, i4)') 'Info: Space group No. ', space_group
1769 write(message(2),'(2a)') 'Info: International: ', trim(symbol)
1770 call messages_info(2, namespace=ions%namespace)
1771
1772 pop_sub(ions_print_spacegroup)
1773 end subroutine ions_print_spacegroup
1774
1775end module ions_oct_m
1776
1777!! Local Variables:
1778!! mode: f90
1779!! coding: utf-8
1780!! 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:1216
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:270
real(real64), parameter, public p_kb
Boltzmann constant in Ha/K.
Definition: global.F90:232
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:683
subroutine ions_update_lattice_vectors(ions, latt, symmetrize)
Regenerate the ions information after update of the lattice vectors.
Definition: ions.F90:1703
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:1001
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:859
real(real64) function, dimension(this%space%dim) ions_global_force(this, time)
Definition: ions.F90:1258
subroutine ions_copy(ions_out, ions_in)
Definition: ions.F90:635
real(real64) function, dimension(this%space%dim) ions_current(this)
Definition: ions.F90:1623
subroutine ions_update_quantity(this, label)
Definition: ions.F90:954
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:1727
subroutine ions_write_xyz(this, fname, append, comment, reduce_coordinates)
Definition: ions.F90:1275
subroutine ions_init_species(ions, factory, print_info)
Definition: ions.F90:523
subroutine ions_finalize(ions)
Definition: ions.F90:1665
subroutine ions_initialize(this)
Definition: ions.F90:940
real(real64) function, dimension(this%space%dim) ions_abs_current(this)
Definition: ions.F90:1644
subroutine ions_read_xyz(this, fname, comment)
Definition: ions.F90:1374
subroutine ions_write_vtk_geometry(this, filename, ascii)
Definition: ions.F90:1519
subroutine ions_rotate(this, from, from2, to)
Definition: ions.F90:1144
subroutine ions_translate(this, xx)
Definition: ions.F90:1128
subroutine ions_symmetrize_atomic_coord(ions)
Symmetrizes atomic coordinates by applying all symmetries.
Definition: ions.F90:1797
subroutine ions_print_spacegroup(ions)
Prints the spacegroup of the system for periodic systems.
Definition: ions.F90:1835
real(real64) function ions_min_distance(this, real_atoms_only)
Definition: ions.F90:1016
subroutine ions_init_random_displacements(ions, T)
create random displacements for positions and velocities
Definition: ions.F90:706
real(real64) function, dimension(this%space%dim) ions_dipole(this, mask)
Definition: ions.F90:1106
subroutine ions_write_bild_forces_file(this, dir, fname)
Definition: ions.F90:1481
subroutine ions_write_crystal(this, dir)
This subroutine creates a crystal by replicating the geometry and writes the result to dir.
Definition: ions.F90:1447
logical function ions_has_time_dependent_species(this)
Definition: ions.F90:1079
subroutine ions_partition(this, mc)
Definition: ions.F90:667
subroutine ions_init_interaction_as_partner(partner, interaction)
Definition: ions.F90:970
subroutine ions_copy_quantities_to_interaction(partner, interaction)
Definition: ions.F90:985
subroutine ions_write_poscar(this, fname, comment)
Writes the positions of the ions in POSCAR format.
Definition: ions.F90:1328
real(real64) function ions_val_charge(this, mask)
Definition: ions.F90:1090
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
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:341
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)