Octopus
ion_dynamics.F90
Go to the documentation of this file.
1!! Copyright (C) 2008 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
23 use iso_c_binding
24 use debug_oct_m
25 use global_oct_m
26 use grid_oct_m
27 use ions_oct_m
28 use, intrinsic :: iso_fortran_env
33 use mpi_oct_m
36 use parser_oct_m
40 use space_oct_m
44 use unit_oct_m
47
48 implicit none
49
50 private
51
52 public :: &
71
72
73 integer, parameter :: &
74 THERMO_NONE = 0, &
75 thermo_scal = 1, &
76 thermo_nh = 2
77
78 type nose_hoover_t
79 private
80 real(real64) :: mass
81 real(real64) :: pos
82 real(real64) :: vel
83 end type nose_hoover_t
84
86 private
87 logical :: move
88 type(tdf_t) :: fx
89 type(tdf_t) :: fy
90 type(tdf_t) :: fz
92
94 private
95 logical :: move_ions
96 logical :: constant_velocity
97 integer :: thermostat
98 real(real64) :: dt
99 real(real64) :: current_temperature
100
101 real(real64), allocatable :: oldforce(:, :)
102
104 real(real64), allocatable :: old_pos(:, :)
105
106 real(real64) :: old_cell_force(3,3)
107 real(real64), public :: stress(3,3)
108 real(real64) :: cell_vel(3,3)
109
111 type(nose_hoover_t) :: nh(1:2)
112 type(tdf_t) :: temperature_function
113
114 logical :: drive_ions
115 type(ion_td_displacement_t), allocatable :: td_displacements(:)
116 type(ions_t), pointer :: ions_t0
117
118 real(real64), public :: ionic_scale
119 end type ion_dynamics_t
120
121 type ion_state_t
122 private
123 real(real64), allocatable :: pos(:, :)
124 real(real64), allocatable :: vel(:, :)
125 real(real64), allocatable :: old_pos(:, :)
126 type(nose_hoover_t) :: nh(1:2)
127 end type ion_state_t
128
129 type cell_state_t
130 private
131 real(real64), allocatable :: pos(:, :)
132 real(real64), allocatable :: vel(:, :)
133 real(real64), allocatable :: old_pos(:, :)
134 end type cell_state_t
135
136contains
137
138 ! ---------------------------------------------------------
139 subroutine ion_dynamics_init(this, namespace, ions)
140 use iso_c_binding, only: c_ptr
141 type(ion_dynamics_t), intent(out) :: this
142 type(namespace_t), intent(in) :: namespace
143 type(ions_t), intent(inout) :: ions
144
145 integer :: i, j, iatom, ierr
146 real(real64) :: xx(ions%space%dim), temperature, sigma, kin1, kin2
147 type(c_ptr) :: random_gen_pointer
148 type(read_coords_info) :: xyz
149 character(len=100) :: temp_function_name
150 logical :: have_velocities
151
152 type(block_t) :: blk
153 integer :: ndisp
154 character(len=200) :: expression
155
156 push_sub(ion_dynamics_init)
157
158 have_velocities = .false.
159 this%drive_ions = .false.
160
161 !%Variable TDIonicTimeScale
162 !%Type float
163 !%Default 1.0
164 !%Section Time-Dependent::Propagation
165 !%Description
166 !% This variable defines the factor between the timescale of ionic
167 !% and electronic movement. It allows reasonably fast
168 !% Born-Oppenheimer molecular-dynamics simulations based on
169 !% Ehrenfest dynamics. The value of this variable is equivalent to
170 !% the role of <math>\mu</math> in Car-Parrinello. Increasing it
171 !% linearly accelerates the time step of the ion
172 !% dynamics, but also increases the deviation of the system from the
173 !% Born-Oppenheimer surface. The default is 1, which means that both
174 !% timescales are the same. Note that a value different than 1
175 !% implies that the electrons will not follow physical behaviour.
176 !%
177 !% According to our tests, values around 10 are reasonable, but it
178 !% will depend on your system, mainly on the width of the gap.
179 !%
180 !% Important: The electronic time step will be the value of
181 !% <tt>TDTimeStep</tt> divided by this variable, so if you have determined an
182 !% optimal electronic time step (that we can call <i>dte</i>), it is
183 !% recommended that you define your time step as:
184 !%
185 !% <tt>TDTimeStep</tt> = <i>dte</i> * <tt>TDIonicTimeScale</tt>
186 !%
187 !% so you will always use the optimal electronic time step
188 !% (<a href=http://arxiv.org/abs/0710.3321>more details</a>).
189 !%End
190 call parse_variable(namespace, 'TDIonicTimeScale', m_one, this%ionic_scale)
192 if (this%ionic_scale <= m_zero) then
193 write(message(1),'(a)') 'Input: TDIonicTimeScale must be positive.'
194 call messages_fatal(1, namespace=namespace)
195 end if
196
197 call messages_print_var_value('TDIonicTimeScale', this%ionic_scale, namespace=namespace)
198
202 !%Variable IonsConstantVelocity
203 !%Type logical
204 !%Default no
205 !%Section Time-Dependent::Propagation
206 !%Description
207 !% (Experimental) If this variable is set to yes, the ions will
208 !% move with a constant velocity given by the initial
209 !% conditions. They will not be affected by any forces.
210 !%End
211 call parse_variable(namespace, 'IonsConstantVelocity', .false., this%constant_velocity)
212 call messages_print_var_value('IonsConstantVelocity', this%constant_velocity, namespace=namespace)
213
214 if (this%constant_velocity) then
215 call messages_experimental('IonsConstantVelocity', namespace=namespace)
216 have_velocities = .true.
217 this%drive_ions = .true.
218 end if
220 !%Variable IonsTimeDependentDisplacements
221 !%Type block
222 !%Section Time-Dependent::Propagation
223 !%Description
224 !% (Experimental) This variable allows you to specify a
225 !% time-dependent function describing the displacement of the ions
226 !% from their equilibrium position: <math>r(t) = r_0 + \Delta
227 !% r(t)</math>. Specify the displacements dx(t), dy(t), dz(t) as
228 !% follows, for some or all of the atoms:
229 !%
230 !% <tt>%IonsTimeDependentDisplacements
231 !% <br>&nbsp;&nbsp; atom_index | "dx(t)" | "dy(t)" | "dz(t)"
232 !% <br>%</tt>
233 !%
234 !% The displacement functions are time-dependent functions and should match one
235 !% of the function names given in the first column of the <tt>TDFunctions</tt> block.
236 !% If this block is set, the ions will not be affected by any forces.
237 !%End
238
239
240 ndisp = 0
241 if (parse_block(namespace, 'IonsTimeDependentDisplacements', blk) == 0) then
242 call messages_experimental("IonsTimeDependentDisplacements", namespace=namespace)
243 ndisp= parse_block_n(blk)
244 safe_allocate(this%td_displacements(1:ions%natoms))
245 this%td_displacements(1:ions%natoms)%move = .false.
246 if (ndisp > 0) this%drive_ions =.true.
247
248 do i = 1, ndisp
249 call parse_block_integer(blk, i-1, 0, iatom)
250 this%td_displacements(iatom)%move = .true.
251
252 call parse_block_string(blk, i-1, 1, expression)
253 call tdf_read(this%td_displacements(iatom)%fx, namespace, trim(expression), ierr)
254 if (ierr /= 0) then
255 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
256 call messages_warning(1, namespace=namespace)
257 end if
258
259
260 call parse_block_string(blk, i-1, 2, expression)
261 call tdf_read(this%td_displacements(iatom)%fy, namespace, trim(expression), ierr)
262 if (ierr /= 0) then
263 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
264 call messages_warning(1, namespace=namespace)
265 end if
266
267 call parse_block_string(blk, i-1, 3, expression)
268 call tdf_read(this%td_displacements(iatom)%fz, namespace, trim(expression), ierr)
269 if (ierr /= 0) then
270 write(message(1),'(3A)') 'Could not find "', trim(expression), '" in the TDFunctions block:'
271 call messages_warning(1, namespace=namespace)
272 end if
273
274 end do
275
276 safe_allocate(this%ions_t0)
277 this%ions_t0 = ions
278
279 end if
280
281
282
283
284 !%Variable Thermostat
285 !%Type integer
286 !%Default none
287 !%Section Time-Dependent::Propagation
288 !%Description
289 !% This variable selects the type of thermostat applied to
290 !% control the ionic temperature.
291 !%Option none 0
292 !% No thermostat is applied. This is the default.
293 !%Option velocity_scaling 1
294 !% Velocities are scaled to control the temperature.
295 !%Option nose_hoover 2
296 !% Nose-Hoover thermostat.
297 !%End
298
299 call parse_variable(namespace, 'Thermostat', thermo_none, this%thermostat)
300 if (.not. varinfo_valid_option('Thermostat', this%thermostat)) call messages_input_error(namespace, 'Thermostat')
301 call messages_print_var_option('Thermostat', this%thermostat, namespace=namespace)
302
303 if (this%thermostat /= thermo_none) then
304
305 have_velocities = .true.
306
307 if (this%drive_ions) then
308 call messages_write('You cannot use a Thermostat and IonsConstantVelocity or IonsTimeDependentDisplacements')
309 call messages_write('at the same time.')
310 call messages_fatal(namespace=namespace)
311 end if
312
313 call messages_experimental('Thermostat', namespace=namespace)
314
315 !%Variable TemperatureFunction
316 !%Type integer
317 !%Default "temperature"
318 !%Section Time-Dependent::Propagation
319 !%Description
320 !% If a thermostat is used, this variable indicates the name of the
321 !% function in the <tt>TDFunctions</tt> block that will be used to control the
322 !% temperature. The values of the temperature are given in
323 !% degrees Kelvin.
324 !%End
325 call parse_variable(namespace, 'TemperatureFunction', 'temperature', temp_function_name)
326
327 call tdf_read(this%temperature_function, namespace, temp_function_name, ierr)
328
329 if (ierr /= 0) then
330 message(1) = "You have enabled a thermostat but Octopus could not find"
331 message(2) = "the '"//trim(temp_function_name)//"' function in the TDFunctions block."
332 call messages_fatal(2, namespace=namespace)
333 end if
334
335 if (this%thermostat == thermo_nh) then
336 !%Variable ThermostatMass
337 !%Type float
338 !%Default 1.0
339 !%Section Time-Dependent::Propagation
340 !%Description
341 !% This variable sets the fictitious mass for the Nose-Hoover
342 !% thermostat.
343 !%End
344 call messages_obsolete_variable(namespace, 'NHMass', 'ThermostatMass')
345
346 call parse_variable(namespace, 'ThermostatMass', m_one, this%nh(1)%mass)
347 this%nh(2)%mass = this%nh(1)%mass
348
349 this%nh(1:2)%pos = m_zero
350 this%nh(1:2)%vel = m_zero
351
352 safe_allocate(this%old_pos(1:ions%space%dim, 1:ions%natoms))
353
354 this%old_pos = ions%pos
355 end if
356
357 end if
358
359 !now initialize velocities
360
361 !%Variable RandomVelocityTemp
362 !%Type float
363 !%Default 0.0
364 !%Section System::Velocities
365 !%Description
366 !% If this variable is present, <tt>Octopus</tt> will assign random
367 !% velocities to the atoms following a Boltzmann distribution with
368 !% temperature given by <tt>RandomVelocityTemp</tt> (in degrees Kelvin).
369 !% The seed for the random number generator can be modified by setting
370 !% <tt>GSL_RNG_SEED</tt> environment variable.
371 !%End
372
373 ! we now load the velocities, either from the temperature, from the input, or from a file
374 if (parse_is_defined(namespace, 'RandomVelocityTemp')) then
375
376 have_velocities = .true.
377
378 if (mpi_grp_is_root(mpi_world)) then
379 call loct_ran_init(random_gen_pointer)
380 call parse_variable(namespace, 'RandomVelocityTemp', m_zero, temperature, unit = unit_kelvin)
381 end if
382
383 do i = 1, ions%natoms
384 !generate the velocities in the root node
385 if (mpi_grp_is_root(mpi_world)) then
386 sigma = sqrt(temperature / ions%mass(i))
387 do j = 1, 3
388 ions%vel(j, i) = loct_ran_gaussian(random_gen_pointer, sigma)
389 end do
390 end if
391 !and send them to the others
392 call mpi_world%bcast(ions%vel(:, i), ions%space%dim, mpi_double_precision, 0)
393 end do
394
395 if (mpi_grp_is_root(mpi_world)) then
396 call loct_ran_end(random_gen_pointer)
397 end if
398
399 call ions%update_kinetic_energy()
400 kin1 = ions%kinetic_energy
401
402 xx = ions%center_of_mass_vel()
403 do i = 1, ions%natoms
404 ions%vel(:, i) = ions%vel(:, i) - xx
405 end do
406
407 call ions%update_kinetic_energy()
408 kin2 = ions%kinetic_energy
409
410 do i = 1, ions%natoms
411 ions%vel(:, i) = sqrt(kin1/kin2)*ions%vel(:, i)
412 end do
413
414 call ions%update_kinetic_energy()
415
416 write(message(1),'(a,f10.4,1x,a)') 'Info: Initial velocities randomly distributed with T =', &
418 write(message(2),'(2x,a,f8.4,1x,a)') '<K> =', &
419 units_from_atomic(units_out%energy, ions%kinetic_energy/ions%natoms), &
420 units_abbrev(units_out%energy)
421 write(message(3),'(2x,a,f8.4,1x,a)') '3/2 k_B T =', &
422 units_from_atomic(units_out%energy, (m_three/m_two)*temperature), &
423 units_abbrev(units_out%energy)
424 call messages_info(3, namespace=namespace)
425
426 else
427 !%Variable XYZVelocities
428 !%Type string
429 !%Section System::Velocities
430 !%Description
431 !% <tt>Octopus</tt> will try to read the starting velocities of the atoms from the XYZ file
432 !% specified by the variable <tt>XYZVelocities</tt>.
433 !% Note that you do not need to specify initial velocities if you are not going
434 !% to perform ion dynamics; if you are going to allow the ions to move but the velocities
435 !% are not specified, they are considered to be null.
436 !% Note: It is important for the velocities to maintain the ordering
437 !% in which the atoms were defined in the coordinates specifications.
438 !%End
439
440 !%Variable XSFVelocities
441 !%Type string
442 !%Section System::Velocities
443 !%Description
444 !% Like <tt>XYZVelocities</tt> but in XCrySDen format, as in <tt>XSFCoordinates</tt>.
445 !%End
446
447 !%Variable PDBVelocities
448 !%Type string
449 !%Section System::Velocities
450 !%Description
451 !% Like <tt>XYZVelocities</tt> but in PDB format, as in <tt>PDBCoordinates</tt>.
452 !%End
453
454 !%Variable Velocities
455 !%Type block
456 !%Section System::Velocities
457 !%Description
458 !% If <tt>XYZVelocities</tt>, <tt>PDBVelocities</tt>, and <tt>XSFVelocities</tt>
459 !% are not present, <tt>Octopus</tt> will try to fetch the initial
460 !% atomic velocities from this block. If this block is not present, <tt>Octopus</tt>
461 !% will set the initial velocities to zero. The format of this block can be
462 !% illustrated by this example:
463 !%
464 !% <tt>%Velocities
465 !% <br>&nbsp;&nbsp;'C' | -1.7 | 0.0 | 0.0
466 !% <br>&nbsp;&nbsp;'O' | &nbsp;1.7 | 0.0 | 0.0
467 !% <br>%</tt>
468 !%
469 !% It describes one carbon and one oxygen moving at the relative
470 !% velocity of 3.4 velocity units.
471 !%
472 !% Note: It is important for the velocities to maintain the ordering
473 !% in which the atoms were defined in the coordinates specifications.
474 !%End
475
476 call read_coords_init(xyz)
477 call read_coords_read('Velocities', xyz, ions%space, namespace)
478 if (xyz%source /= read_coords_err) then
479
480 have_velocities = .true.
481
482 if (ions%natoms /= xyz%n) then
483 write(message(1), '(a,i4,a,i4)') 'I need exactly ', ions%natoms, ' velocities, but I found ', xyz%n
484 call messages_fatal(1, namespace=namespace)
485 end if
486
487 ! copy information and adjust units
488 do i = 1, ions%natoms
489 ions%vel(:, i) = units_to_atomic(units_inp%velocity/units_inp%length, xyz%atom(i)%x(1:ions%space%dim))
490 end do
491
492 call read_coords_end(xyz)
493
494 else
495 ions%vel = m_zero
496 end if
497 end if
498
499 call ions%update_kinetic_energy()
500
501 !%Variable MoveIons
502 !%Type logical
503 !%Section Time-Dependent::Propagation
504 !%Description
505 !% This variable controls whether atoms are moved during a time
506 !% propagation run. The default is yes when the ion velocity is
507 !% set explicitly or implicitly, otherwise is no.
508 !%End
509 call parse_variable(namespace, 'MoveIons', have_velocities, this%move_ions)
510 call messages_print_var_value('MoveIons', this%move_ions, namespace=namespace)
511
512 if (this%move_ions .and. ions%space%periodic_dim == 1) then
513 call messages_input_error(namespace, 'MoveIons', &
514 'Moving ions for a 1D periodic system is not allowed, as forces are incorrect.')
515 end if
516
517 if (ion_dynamics_ions_move(this)) then
518 safe_allocate(this%oldforce(1:ions%space%dim, 1:ions%natoms))
519 end if
520
521 this%old_cell_force = m_zero
522 this%cell_vel = m_zero
523 this%stress = m_zero
524
525
526 pop_sub(ion_dynamics_init)
527 end subroutine ion_dynamics_init
528
529
530 ! ---------------------------------------------------------
531 subroutine ion_dynamics_end(this)
532 type(ion_dynamics_t), intent(inout) :: this
533
534 push_sub(ion_dynamics_end)
535 safe_deallocate_a(this%oldforce)
536
537 if (this%thermostat /= thermo_none) then
538 call tdf_end(this%temperature_function)
539 end if
540
541 if (this%drive_ions .and. allocated(this%td_displacements)) then
542 if (any(this%td_displacements(1:this%ions_t0%natoms)%move)) then
543 ! ions end cannot be called here, otherwise the species are destroyed twice
544 ! call ions_end(this%ions_t0)
545 !AFE_DEALLOCATE_P(this%ions_t0)
546 end if
547 safe_deallocate_a(this%td_displacements)
548 end if
549
550 pop_sub(ion_dynamics_end)
551 end subroutine ion_dynamics_end
552
553
554 ! ---------------------------------------------------------
555 subroutine ion_dynamics_propagate(this, ions, time, dt, namespace)
556 type(ion_dynamics_t), intent(inout) :: this
557 type(ions_t), intent(inout) :: ions
558 real(real64), intent(in) :: time
559 real(real64), intent(in) :: dt
560 type(namespace_t), intent(in) :: namespace
561
562 integer :: iatom!, periodic_dim, idir
563 real(real64) :: dr(3)
564 !real(real64) :: rlattice(3,3)
565
566 if (.not. ion_dynamics_ions_move(this)) return
567
568 push_sub(ion_dynamics_propagate)
569
570 dr = m_zero
571
572 this%dt = dt
573
574
575 ! get the temperature from the tdfunction for the current time
576 if (this%thermostat /= thermo_none) then
577 this%current_temperature = units_to_atomic(unit_kelvin, tdf(this%temperature_function, time))
578
579 if (this%current_temperature < m_zero) then
580 write(message(1), '(a, f10.3, 3a, f10.3, 3a)') &
581 "Negative temperature (", &
582 units_from_atomic(unit_kelvin, this%current_temperature), " ", units_abbrev(unit_kelvin), &
583 ") at time ", &
584 units_from_atomic(units_out%time, time), " ", trim(units_abbrev(units_out%time)), "."
585 call messages_fatal(1, namespace=namespace)
586 end if
587 else
588 this%current_temperature = m_zero
589 end if
590
591 if (this%thermostat /= thermo_nh) then
592 ! integrate using verlet
593 do iatom = 1, ions%natoms
594 if (ions%fixed(iatom)) cycle
595
596 if (.not. this%drive_ions) then
597
598 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom) + &
599 m_half*dt**2 / ions%mass(iatom) * ions%tot_force(:, iatom)
600
601 this%oldforce(:, iatom) = ions%tot_force(:, iatom)
602
603 else if (this%constant_velocity) then
604
605 ions%pos(:, iatom) = ions%pos(:, iatom) + dt*ions%vel(:, iatom)
606
607 else if (allocated(this%td_displacements)) then
608
609 if (this%td_displacements(iatom)%move) then
610
611 dr(1:3)=(/ real(tdf(this%td_displacements(iatom)%fx,time), real64), &
612 real(tdf(this%td_displacements(iatom)%fy,time), real64), &
613 real(tdf(this%td_displacements(iatom)%fz,time), real64) /)
614
615 ions%pos(:, iatom) = this%ions_t0%pos(:, iatom) + dr(1:ions%space%dim)
616
617 end if
618
619 end if
620
621 end do
622
623 else
624 ! for the Nose-Hoover thermostat we use a special integrator
625
626 ! The implementation of the Nose-Hoover thermostat is based on
627 ! Understanding Molecular Simulations by Frenkel and Smit,
628 ! Appendix E, page 540-542.
629
630 call nh_chain(this, ions)
631
632 do iatom = 1, ions%natoms
633 ions%pos(:, iatom) = ions%pos(:, iatom) + m_half*dt*ions%vel(:, iatom)
634 end do
635
636 end if
637
638 ! When the system is periodic in some directions, the atoms might have moved to a an adjacent cell, so we need to move them back to the original cell
639 if (ions%space%is_periodic()) then
640 call ions%fold_atoms_into_cell()
641
642 ! rlattice = ions%latt%rlattice
643
644 ! periodic_dim = ions%space%periodic_dim
645 ! do idir = 1, periodic_dim
646 ! rlattice(1:periodic_dim, idir) = rlattice(1:periodic_dim, idir) &
647 ! + dt*matmul(this%cell_vel(1:periodic_dim, 1:periodic_dim), rlattice(1:periodic_dim, idir)) + &
648 ! M_HALF*dt**2 * matmul(this%stress(1:periodic_dim, 1:periodic_dim), rlattice(1:periodic_dim, idir))
649 ! end do
650 ! this%old_cell_force = this%stress
651 ! call ions%latt%update(rlattice)
652 end if
653
655 end subroutine ion_dynamics_propagate
656
657
658 ! ---------------------------------------------------------
659 subroutine nh_chain(this, ions)
660 type(ion_dynamics_t), intent(inout) :: this
661 type(ions_t), intent(inout) :: ions
662
663 real(real64) :: g1, g2, ss, uk, dt, temp
664
665 push_sub(nh_chain)
666
667 dt = this%dt
668
669 call ions%update_kinetic_energy()
670 uk = ions%kinetic_energy
671
672 temp = this%current_temperature
673
674 g2 = (this%nh(1)%mass*this%nh(1)%vel**2 - temp)/this%nh(2)%mass
675 this%nh(2)%vel = this%nh(2)%vel + g2*dt/m_four
676 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
677
678 g1 = (m_two*uk - m_three*ions%natoms*temp)/this%nh(1)%mass
679 this%nh(1)%vel = this%nh(1)%vel + g1*dt/m_four
680 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
681 this%nh(1)%pos = this%nh(1)%pos + this%nh(1)%vel*dt/m_two
682 this%nh(2)%pos = this%nh(2)%pos + this%nh(2)%vel*dt/m_two
683
684 ss = exp(-this%nh(1)%vel*dt/m_two)
685
686 ions%vel = ss*ions%vel
687
688 uk = uk*ss**2
689
690 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
691 g1 = (m_two*uk - m_three*ions%natoms*temp)/this%nh(1)%mass
692 this%nh(1)%vel = this%nh(1)%vel + g1*dt/m_four
693 this%nh(1)%vel = this%nh(1)%vel*exp(-this%nh(2)%vel*dt/8.0_real64)
694
695 g2 = (this%nh(1)%mass*this%nh(1)%vel**2 - temp)/this%nh(2)%mass
696 this%nh(2)%vel = this%nh(2)%vel + g2*dt/m_four
697
698 pop_sub(nh_chain)
699 end subroutine nh_chain
700
701
702 ! ---------------------------------------------------------
703 subroutine ion_dynamics_propagate_vel(this, ions, atoms_moved)
704 type(ion_dynamics_t), intent(inout) :: this
705 type(ions_t), intent(inout) :: ions
706 logical, optional, intent(out) :: atoms_moved
707
708 integer :: iatom, periodic_dim
709 real(real64) :: scal
710
711 if (.not. ion_dynamics_ions_move(this)) return
712 if (this%drive_ions) return
713
715
716 if (present(atoms_moved)) atoms_moved = this%thermostat == thermo_nh
717
718 if (this%thermostat /= thermo_nh) then
719 ! velocity verlet
720
721 do iatom = 1, ions%natoms
722 if (ions%fixed(iatom)) cycle
723
724 ions%vel(:, iatom) = ions%vel(:, iatom) &
725 + this%dt/ions%mass(iatom) * m_half * (this%oldforce(:, iatom) + &
726 ions%tot_force(:, iatom))
727
728 end do
729
730 else
731 ! the nose-hoover integration
732 do iatom = 1, ions%natoms
733 ions%vel(:, iatom) = ions%vel(:, iatom) + this%dt*ions%tot_force(:, iatom) / ions%mass(iatom)
734 ions%pos(:, iatom) = ions%pos(:, iatom) + m_half*this%dt*ions%vel(:, iatom)
735 end do
736
737 call nh_chain(this, ions)
738
739 end if
740
741 if (this%thermostat == thermo_scal) then
742 scal = sqrt(this%current_temperature/ion_dynamics_temperature(ions))
743
744 ions%vel = scal*ions%vel
745 end if
746
747 if (ions%space%is_periodic()) then
748 periodic_dim = ions%space%periodic_dim
749 this%cell_vel(1:periodic_dim, 1:periodic_dim) = this%cell_vel(1:periodic_dim, 1:periodic_dim) &
750 + this%dt* m_half * (this%old_cell_force(1:periodic_dim, 1:periodic_dim) + this%stress(1:periodic_dim, 1:periodic_dim))
751 end if
754 end subroutine ion_dynamics_propagate_vel
755
756
757 ! ---------------------------------------------------------
759 subroutine ion_dynamics_verlet_step1(ions, q, v, fold, dt)
760 type(ions_t), intent(in) :: ions
761 real(real64), intent(inout) :: q(:, :)
762 real(real64), intent(inout) :: v(:, :)
763 real(real64), intent(in) :: fold(:, :)
764 real(real64), intent(in) :: dt
765
766 integer :: iatom
767
769
770 ! First transform momenta to velocities
771 do iatom = 1, ions%natoms
772 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
773 end do
774
775 ! integrate using verlet
776 do iatom = 1, ions%natoms
777 if (ions%fixed(iatom)) cycle
778 q(iatom, 1:ions%space%dim) = q(iatom, 1:ions%space%dim) + dt * v(iatom, 1:ions%space%dim) + &
779 m_half*dt**2 / ions%mass(iatom) * fold(iatom, 1:ions%space%dim)
780 end do
781
782 ! And back to momenta.
783 do iatom = 1, ions%natoms
784 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
785 end do
786
788 end subroutine ion_dynamics_verlet_step1
789
790
791
792 ! ---------------------------------------------------------
794 subroutine ion_dynamics_verlet_step2(ions, v, fold, fnew, dt)
795 type(ions_t), intent(in) :: ions
796 real(real64), intent(inout) :: v(:, :)
797 real(real64), intent(in) :: fold(:, :)
798 real(real64), intent(in) :: fnew(:, :)
799 real(real64), intent(in) :: dt
800
801 integer :: iatom
802
804
805 ! First transform momenta to velocities
806 do iatom = 1, ions%natoms
807 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) / ions%mass(iatom)
808 end do
809
810 ! velocity verlet
811 do iatom = 1, ions%natoms
812 if (ions%fixed(iatom)) cycle
813 v(iatom, 1:ions%space%dim) = v(iatom, 1:ions%space%dim) &
814 + dt / ions%mass(iatom) * m_half * (fold(iatom, 1:ions%space%dim) + fnew(iatom, 1:ions%space%dim))
815 end do
816
817 ! And back to momenta.
818 do iatom = 1, ions%natoms
819 v(iatom, 1:ions%space%dim) = ions%mass(iatom) * v(iatom, 1:ions%space%dim)
820 end do
821
823 end subroutine ion_dynamics_verlet_step2
824
825
826 ! ---------------------------------------------------------
827 subroutine ion_dynamics_save_state(this, ions, state)
828 type(ion_dynamics_t), intent(in) :: this
829 type(ions_t), intent(in) :: ions
830 type(ion_state_t), intent(out) :: state
831
832 if (.not. ion_dynamics_ions_move(this)) return
833
835
836 safe_allocate(state%pos(1:ions%space%dim, 1:ions%natoms))
837 safe_allocate(state%vel(1:ions%space%dim, 1:ions%natoms))
838
839 state%pos = ions%pos
840 state%vel = ions%vel
841
842 if (this%thermostat == thermo_nh) then
843 safe_allocate(state%old_pos(1:ions%space%dim, 1:ions%natoms))
844 state%old_pos(1:ions%space%dim, 1:ions%natoms) = this%old_pos(1:ions%space%dim, 1:ions%natoms)
845 state%nh(1:2)%pos = this%nh(1:2)%pos
846 state%nh(1:2)%vel = this%nh(1:2)%vel
847 end if
848
850 end subroutine ion_dynamics_save_state
851
853 ! ---------------------------------------------------------
854 subroutine ion_dynamics_restore_state(this, ions, state)
855 type(ion_dynamics_t), intent(inout) :: this
856 type(ions_t), intent(inout) :: ions
857 type(ion_state_t), intent(inout) :: state
858
859 if (.not. ion_dynamics_ions_move(this)) return
860
862
863 ions%pos = state%pos
864 ions%vel = state%vel
865
866 if (this%thermostat == thermo_nh) then
867 this%old_pos(1:ions%space%dim, 1:ions%natoms) = state%old_pos(1:ions%space%dim, 1:ions%natoms)
868 this%nh(1:2)%pos = state%nh(1:2)%pos
869 this%nh(1:2)%vel = state%nh(1:2)%vel
870 safe_deallocate_a(state%old_pos)
871 end if
872
873 safe_deallocate_a(state%pos)
874 safe_deallocate_a(state%vel)
875
877 end subroutine ion_dynamics_restore_state
878
879
880 ! ---------------------------------------------------------
881 logical pure function ion_dynamics_ions_move(this) result(ions_move)
882 type(ion_dynamics_t), intent(in) :: this
883
884 ions_move = this%move_ions
885
886 end function ion_dynamics_ions_move
888
889 ! ---------------------------------------------------------
890 logical pure function ion_dynamics_drive_ions(this) result(drive_ions)
891 type(ion_dynamics_t), intent(in) :: this
892
893 drive_ions = this%drive_ions
894
895 end function ion_dynamics_drive_ions
896
897
898 ! ---------------------------------------------------------
900 real(real64) function ion_dynamics_temperature(ions) result(temperature)
901 type(ions_t), intent(in) :: ions
902
903 integer :: iatom
904 real(real64) :: kinetic_energy
905
906 kinetic_energy = m_zero
907 do iatom = 1, ions%natoms
908 kinetic_energy = kinetic_energy + &
909 m_half * ions%mass(iatom) * sum(ions%vel(:, iatom)**2)
910 end do
911 temperature = m_two/m_three*kinetic_energy/ions%natoms
912
913 end function ion_dynamics_temperature
914
915
916 ! ---------------------------------------------------------
918 logical function ion_dynamics_freeze(this) result(freeze)
919 type(ion_dynamics_t), intent(inout) :: this
920 if (this%move_ions) then
921 this%move_ions = .false.
922 freeze = .true.
923 else
924 freeze = .false.
925 end if
926 end function ion_dynamics_freeze
927
928
929 ! ---------------------------------------------------------
931 subroutine ion_dynamics_unfreeze(this)
932 type(ion_dynamics_t), intent(inout) :: this
933 this%move_ions = .true.
934 end subroutine ion_dynamics_unfreeze
935
936 ! ---------------------------------------------------------
937 subroutine ion_dynamics_dump(this, restart, ierr)
938 type(ion_dynamics_t), intent(in) :: this
939 type(restart_t), intent(in) :: restart
940 integer, intent(out) :: ierr
941
942 push_sub(ion_dynamics_dump)
943
944 if (allocated(this%oldforce)) then
945 call drestart_write_binary(restart, "ion_dynamics_oldforce", size(this%oldforce), &
946 this%oldforce, ierr)
947 end if
948
949 pop_sub(ion_dynamics_dump)
950 end subroutine ion_dynamics_dump
951
952 ! ---------------------------------------------------------
953 subroutine ion_dynamics_load(this, restart, ierr)
954 type(ion_dynamics_t), intent(inout) :: this
955 type(restart_t), intent(in) :: restart
956 integer, intent(out) :: ierr
957
958 push_sub(ion_dynamics_load)
959
960 if (allocated(this%oldforce)) then
961 call drestart_read_binary(restart, "ion_dynamics_oldforce", size(this%oldforce), &
962 this%oldforce, ierr)
963 end if
964
965 pop_sub(ion_dynamics_load)
966 end subroutine ion_dynamics_load
967
968
969 !----------------------------------------------------------
970 subroutine electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
971 type(namespace_t), intent(in) :: namespace
972 type(grid_t), intent(inout) :: gr
973 class(space_t), intent(in) :: space
974 type(poisson_t), intent(inout) :: psolver
975 type(kpoints_t), intent(inout) :: kpoints
976 type(multicomm_t), intent(in) :: mc
977 real(real64), intent(in) :: qtot
978 type(lattice_vectors_t), intent(in) :: new_latt
979
980 integer :: idir
981 real(real64) :: length(1:space%dim)
982
984
985 ! Regenerate the box
986 select type(box => gr%box)
987 type is (box_parallelepiped_t)
988 do idir = 1, space%dim
989 length(idir) = norm2(new_latt%rlattice(1:space%dim, idir))
990 end do
991 call box%regenerate(space%dim, new_latt%rlattice, length, namespace)
992 class default
993 call messages_not_implemented("Grid regeneration for non-parallelepiped boxes", namespace=namespace)
994 end select
995
996
997 call grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
998
999 !Initialize Poisson solvers
1000 call poisson_end(psolver)
1001 call poisson_init(psolver, namespace, space, gr%der, mc, gr%stencil, qtot, verbose=.false.)
1002
1003 call kpoints_lattice_vectors_update(kpoints, new_latt)
1004
1006
1008
1009
1010end module ion_dynamics_oct_m
1012!! Local Variables:
1013!! mode: f90
1014!! coding: utf-8
1015!! End:
Functions to generate random numbers.
Definition: loct_math.F90:329
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
double exp(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_three
Definition: global.F90:190
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public ion_dynamics_verlet_step2(ions, v, fold, fnew, dt)
A bare verlet integrator.
subroutine, public ion_dynamics_dump(this, restart, ierr)
subroutine, public ion_dynamics_verlet_step1(ions, q, v, fold, dt)
A bare verlet integrator.
subroutine nh_chain(this, ions)
subroutine, public ion_dynamics_restore_state(this, ions, state)
subroutine, public ion_dynamics_propagate(this, ions, time, dt, namespace)
subroutine, public ion_dynamics_save_state(this, ions, state)
logical pure function, public ion_dynamics_ions_move(this)
subroutine, public ion_dynamics_unfreeze(this)
Unfreezes the ionic movement.
subroutine, public ion_dynamics_init(this, namespace, ions)
subroutine, public ion_dynamics_propagate_vel(this, ions, atoms_moved)
subroutine, public ion_dynamics_load(this, restart, ierr)
integer, parameter thermo_scal
subroutine, public electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
subroutine, public ion_dynamics_end(this)
integer, parameter thermo_nh
logical function, public ion_dynamics_freeze(this)
Freezes the ionic movement.
logical pure function, public ion_dynamics_drive_ions(this)
real(real64) function, public ion_dynamics_temperature(ions)
This function returns the ionic temperature in energy units.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
integer, parameter, public read_coords_err
for read_coords_info::file_type
subroutine, public read_coords_init(gf)
subroutine, public read_coords_end(gf)
subroutine, public read_coords_read(what, gf, space, namespace)
subroutine, public tdf_end(f)
Definition: tdfunction.F90:974
subroutine, public tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
Definition: tdfunction.F90:218
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_t), public unit_kelvin
For converting energies into temperatures.
type(unit_system_t), public units_inp
the units systems for reading and writing
int true(void)