53#include "../base/base_uses.f90"
59 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'integrator_utils'
63 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: v => null()
64 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: r => null()
65 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: eps => null()
66 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: veps => null()
67 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: h => null()
72 INTEGER,
POINTER :: itimes => null()
73 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: pos => null()
74 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: vel => null()
75 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: shell_pos => null()
76 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: shell_vel => null()
77 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: core_pos => null()
78 REAL(kind=
dp),
POINTER,
DIMENSION(:, :) :: core_vel => null()
79 REAL(kind=
dp) :: max_vel = 0.0_dp, max_vel_sc = 0.0_dp
80 REAL(kind=
dp) :: max_dvel = 0.0_dp, max_dvel_sc = 0.0_dp
81 REAL(kind=
dp) :: max_dr = 0.0_dp, max_dsc = 0.0_dp
82 REAL(kind=
dp) :: arg_r(3) = 0.0_dp, arg_v(3) = 0.0_dp, u(3, 3) = 0.0_dp, e_val(3) = 0.0_dp, s = 0.0_dp, ds = 0.0_dp
83 REAL(kind=
dp),
DIMENSION(3) :: poly_r = 0.0_dp, &
84 poly_v = 0.0_dp, scale_r = 0.0_dp, scale_v = 0.0_dp
88 MODULE PROCEDURE set_particle_set, set_vel
92 MODULE PROCEDURE update_pv_particle_set, update_pv_velocity
96 MODULE PROCEDURE damp_v_particle_set, damp_v_velocity
120 INTEGER :: idim, jdim, natoms
122 natoms =
SIZE(particle_set)
125 cpassert(.NOT.
ASSOCIATED(old))
128 ALLOCATE (old%v(natoms, 3))
130 ALLOCATE (old%r(natoms, 3))
132 ALLOCATE (old%eps(idim, jdim))
134 ALLOCATE (old%veps(idim, jdim))
136 ALLOCATE (old%h(3, 3))
151 IF (
ASSOCIATED(old))
THEN
152 IF (
ASSOCIATED(old%v))
THEN
155 IF (
ASSOCIATED(old%r))
THEN
158 IF (
ASSOCIATED(old%eps))
THEN
161 IF (
ASSOCIATED(old%veps))
THEN
162 DEALLOCATE (old%veps)
164 IF (
ASSOCIATED(old%h))
THEN
184 SUBROUTINE allocate_tmp(md_env, tmp, nparticle, nshell, shell_adiabatic)
187 INTEGER,
INTENT(IN) :: nparticle, nshell
188 LOGICAL,
INTENT(IN) :: shell_adiabatic
190 cpassert(.NOT.
ASSOCIATED(tmp))
196 NULLIFY (tmp%shell_pos)
197 NULLIFY (tmp%shell_vel)
198 NULLIFY (tmp%core_pos)
199 NULLIFY (tmp%core_vel)
203 ALLOCATE (tmp%pos(3, nparticle))
205 ALLOCATE (tmp%vel(3, nparticle))
207 tmp%pos(:, :) = 0.0_dp
208 tmp%vel(:, :) = 0.0_dp
210 IF (shell_adiabatic)
THEN
212 ALLOCATE (tmp%shell_pos(3, nshell))
214 ALLOCATE (tmp%core_pos(3, nshell))
216 ALLOCATE (tmp%shell_vel(3, nshell))
218 ALLOCATE (tmp%core_vel(3, nshell))
220 tmp%shell_pos(:, :) = 0.0_dp
221 tmp%shell_vel(:, :) = 0.0_dp
222 tmp%core_pos(:, :) = 0.0_dp
223 tmp%core_vel(:, :) = 0.0_dp
231 tmp%max_vel_sc = 0.0_dp
232 tmp%max_dvel = 0.0_dp
233 tmp%max_dvel_sc = 0.0_dp
239 CALL get_md_env(md_env=md_env, itimes=tmp%itimes)
259 core_particle_set, para_env, shell_adiabatic, pos, vel, &
263 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, shell_particle_set, &
266 LOGICAL,
INTENT(IN) :: shell_adiabatic
267 LOGICAL,
INTENT(IN),
OPTIONAL :: pos, vel, should_deall_vel
269 LOGICAL :: my_deall, my_pos, my_vel
271 cpassert(
ASSOCIATED(tmp))
275 IF (
PRESENT(pos)) my_pos = pos
276 IF (
PRESENT(vel)) my_vel = vel
277 IF (
PRESENT(should_deall_vel)) my_deall = should_deall_vel
284 IF (shell_adiabatic)
THEN
287 DEALLOCATE (tmp%shell_pos)
288 DEALLOCATE (tmp%core_pos)
295 IF (shell_adiabatic)
THEN
301 IF (shell_adiabatic)
THEN
302 DEALLOCATE (tmp%shell_vel)
303 DEALLOCATE (tmp%core_vel)
305 cpassert(.NOT.
ASSOCIATED(tmp%pos))
306 cpassert(.NOT.
ASSOCIATED(tmp%shell_pos))
307 cpassert(.NOT.
ASSOCIATED(tmp%core_pos))
328 SUBROUTINE get_s_ds(tmp, nparticle_kind, atomic_kind_set, local_particles, particle_set, &
331 INTEGER :: nparticle_kind
337 LOGICAL,
INTENT(IN),
OPTIONAL :: tmpv
339 INTEGER :: iparticle, iparticle_kind, &
340 iparticle_local, nparticle_local
342 REAL(kind=
dp) :: a, b, k, mass, rb
346 IF (
PRESENT(tmpv)) my_tmpv = tmpv
351 DO iparticle_kind = 1, nparticle_kind
352 atomic_kind => atomic_kind_set(iparticle_kind)
354 IF (mass /= 0.0_dp)
THEN
355 nparticle_local = local_particles%n_el(iparticle_kind)
357 DO iparticle_local = 1, nparticle_local
358 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
359 k = k + 0.5_dp*mass*dot_product(tmp%vel(:, iparticle), tmp%vel(:, iparticle))
360 a = a + dot_product(tmp%vel(:, iparticle), particle_set(iparticle)%f(:))
361 b = b + (1.0_dp/mass)*dot_product(particle_set(iparticle)%f(:), particle_set(iparticle)%f(:))
364 DO iparticle_local = 1, nparticle_local
365 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
366 k = k + 0.5_dp*mass*dot_product(particle_set(iparticle)%v(:), particle_set(iparticle)%v(:))
367 a = a + dot_product(particle_set(iparticle)%v(:), particle_set(iparticle)%f(:))
368 b = b + (1.0_dp/mass)*dot_product(particle_set(iparticle)%f(:), particle_set(iparticle)%f(:))
379 tmp%s = (a/b)*(cosh(dt*rb/2.0_dp) - 1) + sinh(dt*rb/2.0_dp)/rb
380 tmp%ds = (a/b)*(sinh(dt*rb/2.0_dp)*rb) + cosh(dt*rb/2.0_dp)
397 SUBROUTINE set_particle_set(old, atomic_kind_set, particle_set, local_particles, cell, npt, char)
404 CHARACTER(LEN=*),
INTENT(IN) :: char
406 INTEGER :: idim, iparticle, iparticle_kind, &
407 iparticle_local, nparticle_kind, &
410 nparticle_kind =
SIZE(atomic_kind_set)
413 DO iparticle_kind = 1, nparticle_kind
414 nparticle_local = local_particles%n_el(iparticle_kind)
415 DO iparticle_local = 1, nparticle_local
416 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
418 old%v(iparticle, idim) = particle_set(iparticle)%v(idim)
419 old%r(iparticle, idim) = particle_set(iparticle)%r(idim)
423 old%eps(:, :) = npt(:, :)%eps
424 old%veps(:, :) = npt(:, :)%v
425 old%h(:, :) = cell%hmat(:, :)
427 DO iparticle_kind = 1, nparticle_kind
428 nparticle_local = local_particles%n_el(iparticle_kind)
429 DO iparticle_local = 1, nparticle_local
430 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
432 particle_set(iparticle)%v(idim) = old%v(iparticle, idim)
433 particle_set(iparticle)%r(idim) = old%r(iparticle, idim)
437 npt(:, :)%eps = old%eps(:, :)
438 npt(:, :)%v = old%veps(:, :)
439 cell%hmat(:, :) = old%h(:, :)
442 END SUBROUTINE set_particle_set
458 SUBROUTINE set_vel(old, atomic_kind_set, particle_set, vel, local_particles, cell, npt, char)
462 REAL(KIND=
dp),
INTENT(INOUT) :: vel(:, :)
466 CHARACTER(LEN=*),
INTENT(IN) :: char
468 INTEGER :: idim, iparticle, iparticle_kind, &
469 iparticle_local, nparticle_kind, &
472 nparticle_kind =
SIZE(atomic_kind_set)
475 DO iparticle_kind = 1, nparticle_kind
476 nparticle_local = local_particles%n_el(iparticle_kind)
477 DO iparticle_local = 1, nparticle_local
478 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
480 old%v(iparticle, idim) = vel(idim, iparticle)
481 old%r(iparticle, idim) = particle_set(iparticle)%r(idim)
485 old%eps(:, :) = npt(:, :)%eps
486 old%veps(:, :) = npt(:, :)%v
487 old%h(:, :) = cell%hmat(:, :)
489 DO iparticle_kind = 1, nparticle_kind
490 nparticle_local = local_particles%n_el(iparticle_kind)
491 DO iparticle_local = 1, nparticle_local
492 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
494 vel(idim, iparticle) = old%v(iparticle, idim)
495 particle_set(iparticle)%r(idim) = old%r(iparticle, idim)
499 npt(:, :)%eps = old%eps(:, :)
500 npt(:, :)%v = old%veps(:, :)
501 cell%hmat(:, :) = old%h(:, :)
504 END SUBROUTINE set_vel
520 SUBROUTINE damp_v_particle_set(molecule_kind_set, molecule_set, &
521 particle_set, local_molecules, gamma1, npt, dt, group)
527 REAL(KIND=
dp),
INTENT(IN) :: gamma1
529 REAL(kind=
dp),
INTENT(IN) :: dt
533 INTEGER :: first_atom, ikind, imol, imol_local, &
534 ipart, last_atom, nmol_local
535 REAL(kind=
dp) :: alpha, ikin, kin, mass, scale
540 DO ikind = 1,
SIZE(molecule_kind_set)
542 nmol_local = local_molecules%n_el(ikind)
543 DO imol_local = 1, nmol_local
544 imol = local_molecules%list(ikind)%array(imol_local)
545 molecule => molecule_set(imol)
548 DO ipart = first_atom, last_atom
549 atomic_kind => particle_set(ipart)%atomic_kind
551 kin = kin + mass*particle_set(ipart)%v(1)* &
552 particle_set(ipart)%v(1)
553 kin = kin + mass*particle_set(ipart)%v(2)* &
554 particle_set(ipart)%v(2)
555 kin = kin + mass*particle_set(ipart)%v(3)* &
556 particle_set(ipart)%v(3)
565 alpha = 2.0_dp*npt%mass*npt%v*npt%v*gamma1*ikin
566 scale = scale*sqrt(1.0_dp + alpha*0.5_dp*dt)
568 DO ikind = 1,
SIZE(molecule_kind_set)
569 nmol_local = local_molecules%n_el(ikind)
570 DO imol_local = 1, nmol_local
571 imol = local_molecules%list(ikind)%array(imol_local)
572 molecule => molecule_set(imol)
575 DO ipart = first_atom, last_atom
576 particle_set(ipart)%v(1) = particle_set(ipart)%v(1)*scale
577 particle_set(ipart)%v(2) = particle_set(ipart)%v(2)*scale
578 particle_set(ipart)%v(3) = particle_set(ipart)%v(3)*scale
583 END SUBROUTINE damp_v_particle_set
600 SUBROUTINE damp_v_velocity(molecule_kind_set, molecule_set, &
601 particle_set, local_molecules, vel, gamma1, npt, dt, group)
607 REAL(KIND=
dp),
INTENT(INOUT) :: vel(:, :)
608 REAL(KIND=
dp),
INTENT(IN) :: gamma1
610 REAL(kind=
dp),
INTENT(IN) :: dt
614 INTEGER :: first_atom, ikind, imol, imol_local, &
615 ipart, last_atom, nmol_local
616 REAL(kind=
dp) :: alpha, ikin, kin, mass, scale
622 DO ikind = 1,
SIZE(molecule_kind_set)
623 nmol_local = local_molecules%n_el(ikind)
624 DO imol_local = 1, nmol_local
625 imol = local_molecules%list(ikind)%array(imol_local)
626 molecule => molecule_set(imol)
629 DO ipart = first_atom, last_atom
630 atomic_kind => particle_set(ipart)%atomic_kind
632 kin = kin + mass*vel(1, ipart)*vel(1, ipart)
633 kin = kin + mass*vel(2, ipart)*vel(2, ipart)
634 kin = kin + mass*vel(3, ipart)*vel(3, ipart)
643 alpha = 2.0_dp*npt%mass*npt%v*npt%v*gamma1*ikin
644 scale = scale*sqrt(1.0_dp + alpha*0.5_dp*dt)
646 DO ikind = 1,
SIZE(molecule_kind_set)
647 nmol_local = local_molecules%n_el(ikind)
648 DO imol_local = 1, nmol_local
649 imol = local_molecules%list(ikind)%array(imol_local)
650 molecule => molecule_set(imol)
653 DO ipart = first_atom, last_atom
654 vel(1, ipart) = vel(1, ipart)*scale
655 vel(2, ipart) = vel(2, ipart)*scale
656 vel(3, ipart) = vel(3, ipart)*scale
661 END SUBROUTINE damp_v_velocity
675 REAL(kind=
dp),
INTENT(IN) :: gamma1, dt
677 REAL(kind=
dp) :: scale
680 scale = scale*exp(-gamma1*0.25_dp*dt)
714 molecule_kind_set, molecule_set, local_molecules, vel, dt, cell, npt, simpar, virial, &
715 vector_v, roll_tol, iroll, infree, first, para_env, u)
725 REAL(kind=
dp),
INTENT(INOUT) :: vel(:, :)
726 REAL(kind=
dp),
INTENT(IN) :: dt
731 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: vector_v
732 REAL(kind=
dp),
INTENT(OUT) :: roll_tol
733 INTEGER,
INTENT(INOUT) :: iroll
734 REAL(kind=
dp),
INTENT(IN) :: infree
735 LOGICAL,
INTENT(INOUT) :: first
737 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: u(:, :)
740 REAL(kind=
dp),
DIMENSION(3, 3) :: pv_kin
744 CALL update_pv(gci, simpar, atomic_kind_set, vel, particle_set, &
745 local_molecules, molecule_set, molecule_kind_set, &
746 local_particles, kin, pv_kin, virial, para_env)
747 CALL update_veps(cell, npt, simpar, pv_kin, kin, virial, infree)
753 SELECT CASE (simpar%ensemble)
755 npt_loc(:, :)%v = 0.0_dp
756 npt_loc(:, :)%mass = 0.0_dp
757 npt_loc(1, 1)%v = npt(1, 1)%v
758 npt_loc(2, 2)%v = npt(1, 1)%v
759 npt_loc(3, 3)%v = npt(1, 1)%v
760 npt_loc(1, 1)%mass = npt(1, 1)%mass
761 npt_loc(2, 2)%mass = npt(1, 1)%mass
762 npt_loc(3, 3)%mass = npt(1, 1)%mass
768 CALL set(old, atomic_kind_set, particle_set, vel, local_particles, cell, npt,
'B')
770 particle_set, vel, dt, simpar, vector_v, npt_loc%v, roll_tol, iroll, &
771 para_env, u, cell, local_particles)
794 SUBROUTINE update_pv_particle_set(gci, simpar, atomic_kind_set, particle_set, &
795 local_molecules, molecule_set, molecule_kind_set, &
796 local_particles, kin, pv_kin, virial, int_group)
805 REAL(KIND=
dp),
INTENT(OUT) :: kin
806 REAL(KIND=
dp),
DIMENSION(3, 3),
INTENT(OUT) :: pv_kin
811 INTEGER :: i, iparticle, iparticle_kind, &
812 iparticle_local, j, nparticle_local
813 REAL(KIND=
dp) :: mass
818 DO iparticle_kind = 1,
SIZE(atomic_kind_set)
819 atomic_kind => atomic_kind_set(iparticle_kind)
821 nparticle_local = local_particles%n_el(iparticle_kind)
822 DO iparticle_local = 1, nparticle_local
823 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
826 pv_kin(i, j) = pv_kin(i, j) + &
827 mass*particle_set(iparticle)%v(i)* &
828 particle_set(iparticle)%v(j)
830 kin = kin + mass*particle_set(iparticle)%v(i)* &
831 particle_set(iparticle)%v(i)
836 CALL int_group%sum(pv_kin)
837 CALL int_group%sum(kin)
840 IF (simpar%constraint)
CALL pv_constraint(gci, local_molecules, &
843 particle_set, virial, &
846 END SUBROUTINE update_pv_particle_set
868 SUBROUTINE update_pv_velocity(gci, simpar, atomic_kind_set, vel, particle_set, &
869 local_molecules, molecule_set, molecule_kind_set, &
870 local_particles, kin, pv_kin, virial, int_group)
874 REAL(KIND=
dp),
INTENT(INOUT) :: vel(:, :)
880 REAL(KIND=
dp),
INTENT(OUT) :: kin
881 REAL(KIND=
dp),
DIMENSION(:, :),
INTENT(OUT) :: pv_kin
886 INTEGER :: i, iparticle, iparticle_kind, &
887 iparticle_local, j, nparticle_local
888 REAL(KIND=
dp) :: mass
893 DO iparticle_kind = 1,
SIZE(atomic_kind_set)
894 atomic_kind => atomic_kind_set(iparticle_kind)
896 nparticle_local = local_particles%n_el(iparticle_kind)
897 DO iparticle_local = 1, nparticle_local
898 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
901 pv_kin(i, j) = pv_kin(i, j) + &
902 mass*vel(i, iparticle)*vel(j, iparticle)
904 kin = kin + mass*vel(i, iparticle)*vel(i, iparticle)
909 CALL int_group%sum(pv_kin)
910 CALL int_group%sum(kin)
913 IF (simpar%constraint)
CALL pv_constraint(gci, local_molecules, &
916 particle_set, virial, &
919 END SUBROUTINE update_pv_velocity
935 SUBROUTINE update_veps(box, npt, simpar, pv_kin, kin, virial, infree, virial_components)
941 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: pv_kin
942 REAL(kind=
dp),
INTENT(IN) :: kin
944 REAL(kind=
dp),
INTENT(IN) :: infree
945 INTEGER,
INTENT(IN),
OPTIONAL :: virial_components
948 REAL(kind=
dp) :: fdotr, trace, v, v0, v0i, vi
949 REAL(kind=
dp),
DIMENSION(3, 3) :: unit
951 cpassert(
ASSOCIATED(box))
966 fdotr = fdotr + virial%pv_virial(ii, ii) + &
967 virial%pv_constraint(ii, ii)
970 npt(:, :)%f = (1.0_dp + (3.0_dp*infree))*kin + fdotr - &
971 3.0_dp*simpar%p_ext*box%deth
974 npt(:, :)%f = virial%pv_virial(:, :) + &
975 pv_kin(:, :) + virial%pv_constraint(:, :) - &
976 unit(:, :)*simpar%p_ext*box%deth + &
977 infree*kin*unit(:, :)
979 trace = npt(1, 1)%f + npt(2, 2)%f + npt(3, 3)%f
981 npt(:, :)%f = trace*unit(:, :)
992 npt(1, 1)%f = virial%pv_virial(1, 1) + &
993 pv_kin(1, 1) + virial%pv_constraint(1, 1) - &
994 simpar%p0*v - simpar%v_shock*simpar%v_shock* &
995 v*v0i*(1._dp - v*v0i) + infree*kin
997 npt(1, 1)%f = virial%pv_virial(1, 1) + &
998 pv_kin(1, 1) + virial%pv_constraint(1, 1) - &
999 simpar%p0*v + infree*kin
1003 npt(1, 1)%f = virial%pv_virial(1, 1) + &
1004 pv_kin(1, 1) + virial%pv_constraint(1, 1) - &
1005 simpar%p0*box%deth + infree*kin
1010 npt(:, :)%v = npt(:, :)%v + &
1011 0.5_dp*simpar%dt*npt(:, :)%f/npt(:, :)%mass
1014 IF (
PRESENT(virial_components))
THEN
1015 SELECT CASE (virial_components)
1019 npt(1, 2)%v = 0.0_dp
1020 npt(1, 3)%v = 0.0_dp
1021 npt(1, 2)%f = 0.0_dp
1022 npt(1, 3)%f = 0.0_dp
1023 npt(2, 1)%v = 0.0_dp
1024 npt(2, 2)%v = 0.0_dp
1025 npt(2, 3)%v = 0.0_dp
1026 npt(2, 1)%f = 0.0_dp
1027 npt(2, 2)%f = 0.0_dp
1028 npt(2, 3)%f = 0.0_dp
1029 npt(3, 1)%v = 0.0_dp
1030 npt(3, 2)%v = 0.0_dp
1031 npt(3, 3)%v = 0.0_dp
1032 npt(3, 1)%f = 0.0_dp
1033 npt(3, 2)%f = 0.0_dp
1034 npt(3, 3)%f = 0.0_dp
1036 npt(1, 1)%v = 0.0_dp
1037 npt(1, 2)%v = 0.0_dp
1038 npt(1, 3)%v = 0.0_dp
1039 npt(1, 1)%f = 0.0_dp
1040 npt(1, 2)%f = 0.0_dp
1041 npt(1, 3)%f = 0.0_dp
1042 npt(2, 1)%v = 0.0_dp
1043 npt(2, 3)%v = 0.0_dp
1044 npt(2, 1)%f = 0.0_dp
1045 npt(2, 3)%f = 0.0_dp
1046 npt(3, 1)%v = 0.0_dp
1047 npt(3, 2)%v = 0.0_dp
1048 npt(3, 3)%v = 0.0_dp
1049 npt(3, 1)%f = 0.0_dp
1050 npt(3, 2)%f = 0.0_dp
1051 npt(3, 3)%f = 0.0_dp
1053 npt(1, 1)%v = 0.0_dp
1054 npt(1, 2)%v = 0.0_dp
1055 npt(1, 3)%v = 0.0_dp
1056 npt(1, 1)%f = 0.0_dp
1057 npt(1, 2)%f = 0.0_dp
1058 npt(1, 3)%f = 0.0_dp
1059 npt(2, 1)%v = 0.0_dp
1060 npt(2, 2)%v = 0.0_dp
1061 npt(2, 3)%v = 0.0_dp
1062 npt(2, 1)%f = 0.0_dp
1063 npt(2, 2)%f = 0.0_dp
1064 npt(2, 3)%f = 0.0_dp
1065 npt(3, 1)%v = 0.0_dp
1066 npt(3, 2)%v = 0.0_dp
1067 npt(3, 1)%f = 0.0_dp
1068 npt(3, 2)%f = 0.0_dp
1070 npt(1, 3)%v = 0.0_dp
1071 npt(1, 3)%f = 0.0_dp
1072 npt(2, 3)%v = 0.0_dp
1073 npt(2, 3)%f = 0.0_dp
1074 npt(3, 1)%v = 0.0_dp
1075 npt(3, 2)%v = 0.0_dp
1076 npt(3, 3)%v = 0.0_dp
1077 npt(3, 1)%f = 0.0_dp
1078 npt(3, 2)%f = 0.0_dp
1079 npt(3, 3)%f = 0.0_dp
1081 npt(1, 2)%v = 0.0_dp
1082 npt(1, 2)%f = 0.0_dp
1083 npt(2, 1)%v = 0.0_dp
1084 npt(2, 2)%v = 0.0_dp
1085 npt(2, 3)%v = 0.0_dp
1086 npt(2, 1)%f = 0.0_dp
1087 npt(2, 2)%f = 0.0_dp
1088 npt(2, 3)%f = 0.0_dp
1089 npt(3, 2)%v = 0.0_dp
1090 npt(3, 2)%f = 0.0_dp
1092 npt(1, 1)%v = 0.0_dp
1093 npt(1, 2)%v = 0.0_dp
1094 npt(1, 3)%v = 0.0_dp
1095 npt(1, 1)%f = 0.0_dp
1096 npt(1, 2)%f = 0.0_dp
1097 npt(1, 3)%f = 0.0_dp
1098 npt(2, 1)%v = 0.0_dp
1099 npt(2, 1)%f = 0.0_dp
1100 npt(3, 1)%v = 0.0_dp
1101 npt(3, 1)%f = 0.0_dp
1125 SUBROUTINE vv_first(tmp, atomic_kind_set, local_particles, particle_set, &
1126 core_particle_set, shell_particle_set, nparticle_kind, &
1127 shell_adiabatic, dt, u, lfixd_list)
1132 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, core_particle_set, &
1134 INTEGER,
INTENT(IN) :: nparticle_kind
1135 LOGICAL,
INTENT(IN) :: shell_adiabatic
1137 REAL(kind=
dp),
DIMENSION(3, 3),
OPTIONAL :: u
1139 OPTIONAL :: lfixd_list
1141 INTEGER :: iparticle, iparticle_kind, &
1142 iparticle_local, nparticle_local, &
1145 REAL(kind=
dp) :: dm, dmc, dms, dsc(3), dvsc(3), &
1146 fac_massc, fac_masss, mass
1150 NULLIFY (atomic_kind, shell)
1151 tmp%max_vel = 0.0_dp
1152 tmp%max_vel_sc = 0.0_dp
1153 tmp%max_dvel = 0.0_dp
1154 tmp%max_dvel_sc = 0.0_dp
1156 tmp%max_dsc = 0.0_dp
1161 DO iparticle_kind = 1, nparticle_kind
1162 atomic_kind => atomic_kind_set(iparticle_kind)
1163 CALL get_atomic_kind(atomic_kind=atomic_kind, mass=mass, shell_active=is_shell)
1164 IF (mass /= 0.0_dp)
THEN
1166 IF (is_shell .AND. shell_adiabatic)
THEN
1168 dms = 0.5_dp*dt/shell%mass_shell
1169 dmc = 0.5_dp*dt/shell%mass_core
1170 fac_masss = shell%mass_shell/mass
1171 fac_massc = shell%mass_core/mass
1172 nparticle_local = local_particles%n_el(iparticle_kind)
1173 IF (
PRESENT(u))
THEN
1174 DO iparticle_local = 1, nparticle_local
1175 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1176 shell_index = particle_set(iparticle)%shell_index
1178 CALL transform_first(shell_particle_set, tmp%shell_pos, tmp%shell_vel, &
1179 shell_index, u, dms, dt, tmp%poly_v, tmp%poly_r, tmp%scale_v, tmp%scale_r)
1182 CALL transform_first(core_particle_set, tmp%core_pos, tmp%core_vel, &
1183 shell_index, u, dmc, dt, tmp%poly_v, tmp%poly_r, tmp%scale_v, tmp%scale_r)
1186 tmp%vel(:, iparticle) = fac_masss*tmp%shell_vel(:, shell_index) + &
1187 fac_massc*tmp%core_vel(:, shell_index)
1188 tmp%pos(:, iparticle) = fac_masss*tmp%shell_pos(:, shell_index) + &
1189 fac_massc*tmp%core_pos(:, shell_index)
1191 tmp%max_vel = max(tmp%max_vel, abs(tmp%vel(1, iparticle)), &
1192 abs(tmp%vel(2, iparticle)), abs(tmp%vel(3, iparticle)))
1193 tmp%max_vel_sc = max(tmp%max_vel_sc, &
1194 abs(tmp%shell_vel(1, shell_index) - tmp%core_vel(1, shell_index)), &
1195 abs(tmp%shell_vel(2, shell_index) - tmp%core_vel(2, shell_index)), &
1196 abs(tmp%shell_vel(3, shell_index) - tmp%core_vel(3, shell_index)))
1197 tmp%max_dr = max(tmp%max_dr, &
1198 abs(particle_set(iparticle)%r(1) - tmp%pos(1, iparticle)), &
1199 abs(particle_set(iparticle)%r(2) - tmp%pos(2, iparticle)), &
1200 abs(particle_set(iparticle)%r(3) - tmp%pos(3, iparticle)))
1201 tmp%max_dvel = max(tmp%max_dvel, &
1202 abs(particle_set(iparticle)%v(1) - tmp%vel(1, iparticle)), &
1203 abs(particle_set(iparticle)%v(2) - tmp%vel(2, iparticle)), &
1204 abs(particle_set(iparticle)%v(3) - tmp%vel(3, iparticle)))
1206 dsc(:) = tmp%shell_pos(:, shell_index) - tmp%core_pos(:, shell_index) - &
1207 shell_particle_set(shell_index)%r(:) + core_particle_set(shell_index)%r(:)
1208 tmp%max_dsc = max(tmp%max_dsc, abs(dsc(1)), abs(dsc(2)), abs(dsc(3)))
1209 dvsc(:) = tmp%shell_vel(:, shell_index) - tmp%core_vel(:, shell_index) - &
1210 shell_particle_set(shell_index)%v(:) + core_particle_set(shell_index)%v(:)
1211 tmp%max_dvel_sc = max(tmp%max_dvel_sc, abs(dvsc(1)), abs(dvsc(2)), abs(dvsc(3)))
1214 DO iparticle_local = 1, nparticle_local
1215 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1216 shell_index = particle_set(iparticle)%shell_index
1217 tmp%shell_vel(:, shell_index) = &
1218 shell_particle_set(shell_index)%v(:)*tmp%scale_v(:)*tmp%scale_v(:) + &
1219 tmp%scale_v(:)*tmp%poly_v(:)*dms*shell_particle_set(shell_index)%f(:)
1220 tmp%shell_pos(:, shell_index) = &
1221 shell_particle_set(shell_index)%r(:)*tmp%scale_r(:)*tmp%scale_r(:) + &
1222 tmp%scale_r(:)*tmp%poly_r(:)*dt*tmp%shell_vel(:, shell_index)
1223 tmp%core_vel(:, shell_index) = &
1224 core_particle_set(shell_index)%v(:)*tmp%scale_v(:)*tmp%scale_v(:) + &
1225 tmp%scale_v(:)*tmp%poly_v(:)*dmc*core_particle_set(shell_index)%f(:)
1226 tmp%core_pos(:, shell_index) = &
1227 core_particle_set(shell_index)%r(:)*tmp%scale_r(:)*tmp%scale_r(:) + &
1228 tmp%scale_r(:)*tmp%poly_r(:)*dt*tmp%core_vel(:, shell_index)
1230 tmp%vel(:, iparticle) = fac_masss*tmp%shell_vel(:, shell_index) + &
1231 fac_massc*tmp%core_vel(:, shell_index)
1232 tmp%pos(:, iparticle) = fac_masss*tmp%shell_pos(:, shell_index) + &
1233 fac_massc*tmp%core_pos(:, shell_index)
1234 tmp%max_vel = max(tmp%max_vel, &
1235 abs(tmp%vel(1, iparticle)), abs(tmp%vel(2, iparticle)), abs(tmp%vel(3, iparticle)))
1236 tmp%max_vel_sc = max(tmp%max_vel_sc, &
1237 abs(tmp%shell_vel(1, shell_index) - tmp%core_vel(1, shell_index)), &
1238 abs(tmp%shell_vel(2, shell_index) - tmp%core_vel(2, shell_index)), &
1239 abs(tmp%shell_vel(3, shell_index) - tmp%core_vel(3, shell_index)))
1240 tmp%max_dr = max(tmp%max_dr, &
1241 abs(particle_set(iparticle)%r(1) - tmp%pos(1, iparticle)), &
1242 abs(particle_set(iparticle)%r(2) - tmp%pos(2, iparticle)), &
1243 abs(particle_set(iparticle)%r(3) - tmp%pos(3, iparticle)))
1244 tmp%max_dvel = max(tmp%max_dvel, &
1245 abs(particle_set(iparticle)%v(1) - tmp%vel(1, iparticle)), &
1246 abs(particle_set(iparticle)%v(2) - tmp%vel(2, iparticle)), &
1247 abs(particle_set(iparticle)%v(3) - tmp%vel(3, iparticle)))
1248 dsc(:) = tmp%shell_pos(:, shell_index) - tmp%core_pos(:, shell_index) - &
1249 shell_particle_set(shell_index)%r(:) + core_particle_set(shell_index)%r(:)
1250 tmp%max_dsc = max(tmp%max_dsc, abs(dsc(1)), abs(dsc(2)), abs(dsc(3)))
1251 dvsc(:) = tmp%shell_vel(:, shell_index) - tmp%core_vel(:, shell_index) - &
1252 shell_particle_set(shell_index)%v(:) + core_particle_set(shell_index)%v(:)
1253 tmp%max_dvel_sc = max(tmp%max_dvel_sc, abs(dvsc(1)), abs(dvsc(2)), abs(dvsc(3)))
1257 nparticle_local = local_particles%n_el(iparticle_kind)
1258 IF (
PRESENT(u))
THEN
1259 DO iparticle_local = 1, nparticle_local
1260 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1262 CALL transform_first(particle_set, tmp%pos, tmp%vel, &
1263 iparticle, u, dm, dt, tmp%poly_v, tmp%poly_r, tmp%scale_v, tmp%scale_r)
1265 tmp%max_vel = max(tmp%max_vel, &
1266 abs(tmp%vel(1, iparticle)), abs(tmp%vel(2, iparticle)), abs(tmp%vel(3, iparticle)))
1267 tmp%max_dr = max(tmp%max_dr, abs(particle_set(iparticle)%r(1) - tmp%pos(1, iparticle)), &
1268 abs(particle_set(iparticle)%r(2) - tmp%pos(2, iparticle)), &
1269 abs(particle_set(iparticle)%r(3) - tmp%pos(3, iparticle)))
1270 tmp%max_dvel = max(tmp%max_dvel, &
1271 abs(particle_set(iparticle)%v(1) - tmp%vel(1, iparticle)), &
1272 abs(particle_set(iparticle)%v(2) - tmp%vel(2, iparticle)), &
1273 abs(particle_set(iparticle)%v(3) - tmp%vel(3, iparticle)))
1276 DO iparticle_local = 1, nparticle_local
1277 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1278 tmp%vel(1, iparticle) = &
1279 particle_set(iparticle)%v(1)*tmp%scale_v(1)*tmp%scale_v(1) + &
1280 tmp%scale_v(1)*tmp%poly_v(1)*dm*particle_set(iparticle)%f(1)
1281 tmp%vel(2, iparticle) = &
1282 particle_set(iparticle)%v(2)*tmp%scale_v(2)*tmp%scale_v(2) + &
1283 tmp%scale_v(2)*tmp%poly_v(2)*dm*particle_set(iparticle)%f(2)
1284 tmp%vel(3, iparticle) = &
1285 particle_set(iparticle)%v(3)*tmp%scale_v(3)*tmp%scale_v(3) + &
1286 tmp%scale_v(3)*tmp%poly_v(3)*dm*particle_set(iparticle)%f(3)
1287 tmp%pos(1, iparticle) = &
1288 particle_set(iparticle)%r(1)*tmp%scale_r(1)*tmp%scale_r(1) + &
1289 tmp%scale_r(1)*tmp%poly_r(1)*dt*tmp%vel(1, iparticle)
1290 tmp%pos(2, iparticle) = &
1291 particle_set(iparticle)%r(2)*tmp%scale_r(2)*tmp%scale_r(2) + &
1292 tmp%scale_r(2)*tmp%poly_r(2)*dt*tmp%vel(2, iparticle)
1293 tmp%pos(3, iparticle) = &
1294 particle_set(iparticle)%r(3)*tmp%scale_r(3)*tmp%scale_r(3) + &
1295 tmp%scale_r(3)*tmp%poly_r(3)*dt*tmp%vel(3, iparticle)
1298 IF (
PRESENT(lfixd_list))
THEN
1299 IF (any(lfixd_list(:)%ifixd_index == iparticle))
THEN
1300 tmp%pos(1, iparticle) = particle_set(iparticle)%r(1)
1301 tmp%pos(2, iparticle) = particle_set(iparticle)%r(2)
1302 tmp%pos(3, iparticle) = particle_set(iparticle)%r(3)
1306 tmp%max_vel = max(tmp%max_vel, &
1307 abs(tmp%vel(1, iparticle)), abs(tmp%vel(2, iparticle)), abs(tmp%vel(3, iparticle)))
1308 tmp%max_dr = max(tmp%max_dr, &
1309 abs(particle_set(iparticle)%r(1) - tmp%pos(1, iparticle)), &
1310 abs(particle_set(iparticle)%r(2) - tmp%pos(2, iparticle)), &
1311 abs(particle_set(iparticle)%r(3) - tmp%pos(3, iparticle)))
1312 tmp%max_dvel = max(tmp%max_dvel, &
1313 abs(particle_set(iparticle)%v(1) - tmp%vel(1, iparticle)), &
1314 abs(particle_set(iparticle)%v(2) - tmp%vel(2, iparticle)), &
1315 abs(particle_set(iparticle)%v(3) - tmp%vel(3, iparticle)))
1340 SUBROUTINE vv_second(tmp, atomic_kind_set, local_particles, particle_set, &
1341 core_particle_set, shell_particle_set, nparticle_kind, &
1342 shell_adiabatic, dt, u)
1347 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, core_particle_set, &
1349 INTEGER,
INTENT(IN) :: nparticle_kind
1350 LOGICAL,
INTENT(IN) :: shell_adiabatic
1352 REAL(kind=
dp),
DIMENSION(3, 3),
OPTIONAL :: u
1354 INTEGER :: iparticle, iparticle_kind, &
1355 iparticle_local, nparticle_local, &
1358 REAL(kind=
dp) :: dm, dmc, dms, fac_massc, fac_masss, mass
1362 NULLIFY (atomic_kind, shell)
1365 DO iparticle_kind = 1, nparticle_kind
1366 atomic_kind => atomic_kind_set(iparticle_kind)
1368 shell_active=is_shell)
1369 IF (mass /= 0.0_dp)
THEN
1371 IF (is_shell .AND. shell_adiabatic)
THEN
1373 dms = 0.5_dp*dt/shell%mass_shell
1374 dmc = 0.5_dp*dt/shell%mass_core
1375 fac_masss = shell%mass_shell/mass
1376 fac_massc = shell%mass_core/mass
1377 nparticle_local = local_particles%n_el(iparticle_kind)
1378 IF (
PRESENT(u))
THEN
1379 DO iparticle_local = 1, nparticle_local
1380 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1381 shell_index = particle_set(iparticle)%shell_index
1383 CALL transform_second(shell_particle_set, tmp%shell_vel, shell_index, &
1384 u, dms, tmp%poly_v, tmp%scale_v)
1387 CALL transform_second(core_particle_set, tmp%core_vel, shell_index, &
1388 u, dmc, tmp%poly_v, tmp%scale_v)
1391 tmp%vel(1, iparticle) = fac_masss*tmp%shell_vel(1, shell_index) + &
1392 fac_massc*tmp%core_vel(1, shell_index)
1393 tmp%vel(2, iparticle) = fac_masss*tmp%shell_vel(2, shell_index) + &
1394 fac_massc*tmp%core_vel(2, shell_index)
1395 tmp%vel(3, iparticle) = fac_masss*tmp%shell_vel(3, shell_index) + &
1396 fac_massc*tmp%core_vel(3, shell_index)
1399 DO iparticle_local = 1, nparticle_local
1400 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1401 shell_index = particle_set(iparticle)%shell_index
1402 tmp%shell_vel(1, shell_index) = &
1403 tmp%shell_vel(1, shell_index)*tmp%scale_v(1)*tmp%scale_v(1) + &
1404 tmp%scale_v(1)*tmp%poly_v(1)*dms*shell_particle_set(shell_index)%f(1)
1405 tmp%shell_vel(2, shell_index) = &
1406 tmp%shell_vel(2, shell_index)*tmp%scale_v(2)*tmp%scale_v(2) + &
1407 tmp%scale_v(2)*tmp%poly_v(2)*dms*shell_particle_set(shell_index)%f(2)
1408 tmp%shell_vel(3, shell_index) = &
1409 tmp%shell_vel(3, shell_index)*tmp%scale_v(3)*tmp%scale_v(3) + &
1410 tmp%scale_v(3)*tmp%poly_v(3)*dms*shell_particle_set(shell_index)%f(3)
1411 tmp%core_vel(1, shell_index) = &
1412 tmp%core_vel(1, shell_index)*tmp%scale_v(1)*tmp%scale_v(1) + &
1413 tmp%scale_v(1)*tmp%poly_v(1)*dmc*core_particle_set(shell_index)%f(1)
1414 tmp%core_vel(2, shell_index) = &
1415 tmp%core_vel(2, shell_index)*tmp%scale_v(2)*tmp%scale_v(2) + &
1416 tmp%scale_v(2)*tmp%poly_v(2)*dmc*core_particle_set(shell_index)%f(2)
1417 tmp%core_vel(3, shell_index) = &
1418 tmp%core_vel(3, shell_index)*tmp%scale_v(3)*tmp%scale_v(3) + &
1419 tmp%scale_v(3)*tmp%poly_v(3)*dmc*core_particle_set(shell_index)%f(3)
1421 tmp%vel(1, iparticle) = fac_masss*tmp%shell_vel(1, shell_index) + &
1422 fac_massc*tmp%core_vel(1, shell_index)
1423 tmp%vel(2, iparticle) = fac_masss*tmp%shell_vel(2, shell_index) + &
1424 fac_massc*tmp%core_vel(2, shell_index)
1425 tmp%vel(3, iparticle) = fac_masss*tmp%shell_vel(3, shell_index) + &
1426 fac_massc*tmp%core_vel(3, shell_index)
1430 nparticle_local = local_particles%n_el(iparticle_kind)
1431 IF (
PRESENT(u))
THEN
1432 DO iparticle_local = 1, nparticle_local
1433 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1434 CALL transform_second(particle_set, tmp%vel, iparticle, &
1435 u, dm, tmp%poly_v, tmp%scale_v)
1438 DO iparticle_local = 1, nparticle_local
1439 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
1440 tmp%vel(1, iparticle) = &
1441 tmp%vel(1, iparticle)*tmp%scale_v(1)*tmp%scale_v(1) + &
1442 tmp%scale_v(1)*tmp%poly_v(1)*dm*particle_set(iparticle)%f(1)
1443 tmp%vel(2, iparticle) = &
1444 tmp%vel(2, iparticle)*tmp%scale_v(2)*tmp%scale_v(2) + &
1445 tmp%scale_v(2)*tmp%poly_v(2)*dm*particle_set(iparticle)%f(2)
1446 tmp%vel(3, iparticle) = &
1447 tmp%vel(3, iparticle)*tmp%scale_v(3)*tmp%scale_v(3) + &
1448 tmp%scale_v(3)*tmp%poly_v(3)*dm*particle_set(iparticle)%f(3)
1475 SUBROUTINE transform_first(particle_set, pos, vel, index, u, dm, dt, poly_v, &
1476 poly_r, scale_v, scale_r)
1479 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: pos, vel
1480 INTEGER,
INTENT(IN) :: index
1481 REAL(kind=
dp),
INTENT(IN) :: u(3, 3), dm, dt, poly_v(3), poly_r(3), &
1482 scale_v(3), scale_r(3)
1484 REAL(kind=
dp),
DIMENSION(3) :: uf, ur, uv
1486 ur = matmul(transpose(u), particle_set(index)%r(:))
1487 uv = matmul(transpose(u), particle_set(index)%v(:))
1488 uf = matmul(transpose(u), particle_set(index)%f(:))
1490 uv(1) = uv(1)*scale_v(1)*scale_v(1) + uf(1)*scale_v(1)*poly_v(1)*dm
1491 uv(2) = uv(2)*scale_v(2)*scale_v(2) + uf(2)*scale_v(2)*poly_v(2)*dm
1492 uv(3) = uv(3)*scale_v(3)*scale_v(3) + uf(3)*scale_v(3)*poly_v(3)*dm
1494 ur(1) = ur(1)*scale_r(1)*scale_r(1) + uv(1)*scale_r(1)*poly_r(1)*dt
1495 ur(2) = ur(2)*scale_r(2)*scale_r(2) + uv(2)*scale_r(2)*poly_r(2)*dt
1496 ur(3) = ur(3)*scale_r(3)*scale_r(3) + uv(3)*scale_r(3)*poly_r(3)*dt
1498 pos(:, index) = matmul(u, ur)
1499 vel(:, index) = matmul(u, uv)
1501 END SUBROUTINE transform_first
1517 SUBROUTINE transform_second(particle_set, vel, index, u, dm, poly_v, scale_v)
1520 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: vel
1521 INTEGER,
INTENT(IN) :: index
1522 REAL(kind=
dp),
INTENT(IN) :: u(3, 3), dm, poly_v(3), scale_v(3)
1524 REAL(kind=
dp),
DIMENSION(3) :: uf, uv
1526 uv = matmul(transpose(u), vel(:, index))
1527 uf = matmul(transpose(u), particle_set(index)%f(:))
1529 uv(1) = uv(1)*scale_v(1)*scale_v(1) + uf(1)*scale_v(1)*poly_v(1)*dm
1530 uv(2) = uv(2)*scale_v(2)*scale_v(2) + uf(2)*scale_v(2)*poly_v(2)*dm
1531 uv(3) = uv(3)*scale_v(3)*scale_v(3) + uf(3)*scale_v(3)*poly_v(3)*dm
1533 vel(:, index) = matmul(u, uv)
1535 END SUBROUTINE transform_second
1557 particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, npt)
1561 REAL(kind=
dp),
INTENT(INOUT) :: dt
1566 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, core_particle_set, &
1568 INTEGER,
INTENT(IN) :: nparticle_kind
1569 LOGICAL,
INTENT(IN) :: shell_adiabatic
1572 INTEGER,
SAVE :: itime = 0
1573 REAL(kind=
dp) :: dt_fact, dt_fact_f, dt_fact_old, &
1574 dt_fact_sc, dt_fact_v, dt_old
1575 REAL(kind=
dp),
POINTER :: time
1583 dt_fact_old = simpar%dt_fact
1584 simpar%dt_fact = 1.0_dp
1585 NULLIFY (thermostats)
1588 CALL para_env%max(tmp%max_dr)
1589 IF (tmp%max_dr > simpar%dr_tol)
THEN
1590 CALL para_env%max(tmp%max_dvel)
1591 dt_fact = sqrt(simpar%dr_tol*dt/tmp%max_dvel)/dt
1594 IF (shell_adiabatic)
THEN
1595 CALL para_env%max(tmp%max_dsc)
1596 IF (tmp%max_dsc > simpar%dsc_tol)
THEN
1597 CALL para_env%max(tmp%max_dvel_sc)
1598 dt_fact_sc = sqrt(simpar%dsc_tol*dt/tmp%max_dvel_sc)/dt
1602 dt_fact_f = min(dt_fact, dt_fact_sc, 1.0_dp)
1603 IF (dt_fact_f < 1.0_dp)
THEN
1604 IF (dt_fact_f < 0.1_dp) dt_fact_f = 0.1_dp
1607 simpar%dt_fact = dt_fact_f
1608 CALL rescaled_vv_first(tmp, dt, simpar, atomic_kind_set, local_particles, &
1609 particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, npt)
1614 CALL para_env%max(tmp%max_dr)
1615 IF (tmp%max_dr > simpar%dr_tol)
THEN
1616 CALL para_env%max(tmp%max_vel)
1617 dt_fact = simpar%dr_tol/tmp%max_vel/dt
1620 IF (shell_adiabatic)
THEN
1621 CALL para_env%max(tmp%max_dsc)
1622 IF (tmp%max_dsc > simpar%dsc_tol)
THEN
1623 CALL para_env%max(tmp%max_vel_sc)
1624 dt_fact_sc = simpar%dsc_tol/tmp%max_vel_sc/dt
1627 dt_fact_v = min(dt_fact, dt_fact_sc, 1.0_dp)
1629 IF (dt_fact_v < 1.0_dp)
THEN
1630 IF (dt_fact_v < 0.1_dp) dt_fact_v = 0.1_dp
1633 simpar%dt_fact = dt_fact_f*dt_fact_v
1634 CALL rescaled_vv_first(tmp, dt, simpar, atomic_kind_set, local_particles, &
1635 particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, npt)
1638 simpar%dt_fact = dt_fact_f*dt_fact_v
1639 IF (simpar%dt_fact < 1.0_dp)
THEN
1640 CALL get_md_env(md_env, t=time, thermostats=thermostats)
1641 time = time - dt_old + dt_old*simpar%dt_fact
1642 IF (
ASSOCIATED(thermostats))
THEN
1667 SUBROUTINE rescaled_vv_first(tmp, dt, simpar, atomic_kind_set, local_particles, &
1668 particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, npt)
1671 REAL(kind=
dp),
INTENT(IN) :: dt
1675 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, core_particle_set, &
1677 INTEGER,
INTENT(IN) :: nparticle_kind
1678 LOGICAL,
INTENT(IN) :: shell_adiabatic
1681 REAL(kind=
dp),
PARAMETER :: e2 = 1.0_dp/6.0_dp, e4 = e2/20.0_dp, &
1682 e6 = e4/42.0_dp, e8 = e6/72.0_dp
1684 REAL(kind=
dp) :: arg_r(3), arg_v(3), e_val(3), infree, &
1687 infree = 1.0_dp/real(simpar%nfree,
dp)
1693 SELECT CASE (simpar%ensemble)
1695 CALL vv_first(tmp, atomic_kind_set, local_particles, particle_set, &
1696 core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, dt)
1698 tmp%scale_v(1:3) = sqrt(1.0_dp/tmp%ds)
1699 tmp%poly_v(1:3) = 2.0_dp*tmp%s/sqrt(tmp%ds)/dt
1700 CALL vv_first(tmp, atomic_kind_set, local_particles, particle_set, &
1701 core_particle_set, shell_particle_set, nparticle_kind, &
1702 shell_adiabatic, dt)
1705 arg_r = arg_r*simpar%dt_fact*simpar%dt_fact
1706 tmp%poly_r(1:3) = 1.0_dp + e2*arg_r(1) + e4*arg_r(1)*arg_r(1) + e6*arg_r(1)**3 + e8*arg_r(1)**4
1707 arg_v = arg_v*simpar%dt_fact*simpar%dt_fact
1708 tmp%poly_v(1:3) = 1.0_dp + e2*arg_v(1) + e4*arg_v(1)*arg_v(1) + e6*arg_v(1)**3 + e8*arg_v(1)**4
1709 tmp%scale_r(1:3) = exp(0.5_dp*dt*npt(1, 1)%v)
1710 tmp%scale_v(1:3) = exp(-0.25_dp*dt*npt(1, 1)%v* &
1711 (1.0_dp + 3.0_dp*infree))
1712 CALL vv_first(tmp, atomic_kind_set, local_particles, particle_set, &
1713 core_particle_set, shell_particle_set, nparticle_kind, &
1714 shell_adiabatic, dt)
1717 trvg = npt(1, 1)%v + npt(2, 2)%v + npt(3, 3)%v
1718 arg_r(:) = arg_r(:)*simpar%dt_fact*simpar%dt_fact
1719 tmp%poly_r = 1._dp + e2*arg_r + e4*arg_r*arg_r + e6*arg_r**3 + e8*arg_r**4
1720 tmp%scale_r(:) = exp(0.5_dp*dt*e_val(:))
1721 arg_v(:) = arg_v(:)*simpar%dt_fact*simpar%dt_fact
1722 tmp%poly_v = 1.0_dp + e2*arg_v + e4*arg_v*arg_v + e6*arg_v**3 + e8*arg_v**4
1723 tmp%scale_v(:) = exp(-0.25_dp*dt*( &
1724 e_val(:) + trvg*infree))
1726 CALL vv_first(tmp, atomic_kind_set, local_particles, particle_set, &
1727 core_particle_set, shell_particle_set, nparticle_kind, &
1728 shell_adiabatic, dt, u)
1731 arg_r = arg_r*simpar%dt_fact*simpar%dt_fact
1732 tmp%poly_r(1) = 1._dp + e2*arg_r(1) + e4*arg_r(1)*arg_r(1) + e6*arg_r(1)**3 + e8*arg_r(1)**4
1733 arg_v(2) = arg_v(2)*simpar%dt_fact*simpar%dt_fact
1734 arg_v(1) = arg_v(1)*simpar%dt_fact*simpar%dt_fact
1735 tmp%poly_v(1) = 1._dp + e2*arg_v(1) + e4*arg_v(1)*arg_v(1) + e6*arg_v(1)**3 + e8*arg_v(1)**4
1736 tmp%poly_v(2) = 1._dp + e2*arg_v(2) + e4*arg_v(2)*arg_v(2) + e6*arg_v(2)**3 + e8*arg_v(2)**4
1737 tmp%poly_v(3) = 1._dp + e2*arg_v(2) + e4*arg_v(2)*arg_v(2) + e6*arg_v(2)**3 + e8*arg_v(2)**4
1738 tmp%scale_r(1) = exp(0.5_dp*dt*npt(1, 1)%v)
1739 tmp%scale_v(1) = exp(-0.25_dp*dt*npt(1, 1)%v* &
1741 tmp%scale_v(2) = exp(-0.25_dp*dt*npt(1, 1)%v*infree)
1742 tmp%scale_v(3) = exp(-0.25_dp*dt*npt(1, 1)%v*infree)
1743 CALL vv_first(tmp, atomic_kind_set, local_particles, particle_set, &
1744 core_particle_set, shell_particle_set, nparticle_kind, &
1745 shell_adiabatic, dt)
1749 END SUBROUTINE rescaled_vv_first
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Barostat structure: module containing barostat available for MD.
integer, parameter, public do_clv_y
integer, parameter, public do_clv_xyz
integer, parameter, public do_clv_yz
integer, parameter, public do_clv_xy
integer, parameter, public do_clv_z
integer, parameter, public do_clv_xz
integer, parameter, public do_clv_x
Handles all functions related to the CELL.
Contains routines useful for the application of constraints during MD.
subroutine, public pv_constraint(gci, local_molecules, molecule_set, molecule_kind_set, particle_set, virial, group)
...
subroutine, public rattle_roll_control(gci, local_molecules, molecule_set, molecule_kind_set, particle_set, vel, dt, simpar, vector, veps, roll_tol, iroll, para_env, u, cell, local_particles)
...
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
Lumps all possible extended system variables into one type for easy access and passing.
logical, parameter, public debug_uniaxial_limit
logical, parameter, public debug_isotropic_limit
Provides integrator utility routines for the integrators.
subroutine, public variable_timestep(md_env, tmp, dt, simpar, para_env, atomic_kind_set, local_particles, particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, npt)
Compute the timestep rescaling factor.
subroutine, public rattle_roll_setup(old, gci, atomic_kind_set, particle_set, local_particles, molecule_kind_set, molecule_set, local_molecules, vel, dt, cell, npt, simpar, virial, vector_v, roll_tol, iroll, infree, first, para_env, u)
update veps using multiplier obtained from SHAKE
subroutine, public allocate_tmp(md_env, tmp, nparticle, nshell, shell_adiabatic)
allocate temporary variables to store positions and velocities used by the velocity-verlet integrator
subroutine, public vv_second(tmp, atomic_kind_set, local_particles, particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, dt, u)
Second half of the velocity-verlet algorithm : update velocity by half step using the new forces.
subroutine, public allocate_old(old, particle_set, npt)
...
subroutine, public update_dealloc_tmp(tmp, particle_set, shell_particle_set, core_particle_set, para_env, shell_adiabatic, pos, vel, should_deall_vel)
update positions and deallocate temporary variable
elemental subroutine, public damp_veps(npt, gamma1, dt)
provides damping for barostat via nph_uniaxial_damped dynamics
subroutine, public update_veps(box, npt, simpar, pv_kin, kin, virial, infree, virial_components)
Routine to compute veps.
subroutine, public get_s_ds(tmp, nparticle_kind, atomic_kind_set, local_particles, particle_set, dt, para_env, tmpv)
...
subroutine, public vv_first(tmp, atomic_kind_set, local_particles, particle_set, core_particle_set, shell_particle_set, nparticle_kind, shell_adiabatic, dt, u, lfixd_list)
First half of the velocity-verlet algorithm : update velocity by half step and positions by full step...
subroutine, public deallocate_old(old)
...
Defines the basic variable types.
integer, parameter, public dp
subroutine, public get_md_env(md_env, itimes, constant, used_time, cell, simpar, npt, force_env, para_env, reftraj, t, init, first_time, fe_env, thermostats, barostat, thermostat_coeff, thermostat_part, thermostat_shell, thermostat_baro, thermostat_fast, thermostat_slow, md_ener, averages, thermal_regions, ehrenfest_md)
get components of MD environment type
Interface to the message passing library MPI.
Define the molecule kind structure types and the corresponding functionality.
Define the data structure for the molecule information.
subroutine, public get_molecule(molecule, molecule_kind, lmi, lci, lg3x3, lg4x6, lcolv, first_atom, last_atom, first_shell, last_shell)
Get components from a molecule data set.
Define the data structure for the particle information.
subroutine, public update_particle_set(particle_set, int_group, pos, vel, for, add)
...
Type for storing MD parameters.
Thermostat structure: module containing thermostat available for MD.
subroutine, public set_thermostats(thermostats, dt_fact)
access internal structures of thermostats
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
structure to store local (to a processor) ordered lists of integers.
stores all the informations relevant to an mpi environment
Simulation parameter type for molecular dynamics.