1 |
chuckv |
4 |
module dynamics_nose_hoover |
2 |
|
|
use definitions, ONLY : DP,NDIM |
3 |
|
|
use constants |
4 |
|
|
use simulation |
5 |
|
|
use parameter |
6 |
|
|
use force_module, ONLY : calc_forces |
7 |
|
|
use force_utilities, ONLY : check |
8 |
|
|
use velocity, ONLY : scale_velocities |
9 |
|
|
use dynamics_utilities, ONLY : evolve_time, sim_status, init_dynamics_loop, & |
10 |
|
|
DT2,DTSQRT,init_status |
11 |
|
|
#ifdef MPI |
12 |
|
|
use mpi_module, ONLY : mpi_allreduce,mpi_comm_world,mpi_err,& |
13 |
|
|
mpi_double_precision, mpi_sum |
14 |
|
|
#endif |
15 |
|
|
|
16 |
|
|
implicit none |
17 |
|
|
PRIVATE |
18 |
|
|
|
19 |
|
|
|
20 |
|
|
|
21 |
|
|
public :: nose_hoover_dynamics |
22 |
|
|
|
23 |
|
|
type (sim_status) :: nose_status |
24 |
|
|
|
25 |
|
|
|
26 |
|
|
contains |
27 |
|
|
|
28 |
|
|
|
29 |
|
|
subroutine nose_hoover_dynamics() |
30 |
|
|
use velocity, ONLY : calc_temp |
31 |
|
|
|
32 |
|
|
|
33 |
|
|
|
34 |
|
|
logical :: update_nlist |
35 |
|
|
logical :: do_pot |
36 |
|
|
logical :: not_done |
37 |
|
|
logical :: nmflag |
38 |
|
|
|
39 |
|
|
integer :: step, i |
40 |
|
|
|
41 |
|
|
call init_status(nose_status) |
42 |
|
|
call init_dynamics_loop( nose_status ) |
43 |
|
|
|
44 |
|
|
! Start the main simulation loop. |
45 |
|
|
|
46 |
|
|
do_pot = .true. |
47 |
|
|
not_done = .true. |
48 |
|
|
nmflag = .false. |
49 |
|
|
step = 0 |
50 |
|
|
|
51 |
|
|
dynamics: do |
52 |
|
|
if (.not.not_done) exit dynamics |
53 |
|
|
|
54 |
|
|
step = step + 1 |
55 |
|
|
|
56 |
|
|
|
57 |
|
|
if (checktemp.and.(mod(step,check_temp_steps).eq.0)) then |
58 |
|
|
call calc_temp(nose_status%temp) |
59 |
|
|
if (dabs(nose_status%temp-target_temp).gt.therm_variance) then |
60 |
|
|
call scale_velocities(target_temp) |
61 |
|
|
end if |
62 |
|
|
end if |
63 |
|
|
|
64 |
|
|
call hooverb() |
65 |
|
|
|
66 |
|
|
call check(update_nlist) |
67 |
|
|
|
68 |
|
|
|
69 |
|
|
if (do_pot) then |
70 |
|
|
call calc_forces(update_nlist,nmflag,pe = nose_status%pot_e) |
71 |
|
|
else |
72 |
|
|
call calc_forces(update_nlist,nmflag) |
73 |
|
|
endif |
74 |
|
|
|
75 |
|
|
|
76 |
|
|
call hoovera() |
77 |
|
|
|
78 |
|
|
call evolve_time(step,nose_status,do_pot,not_done) |
79 |
|
|
|
80 |
|
|
|
81 |
|
|
#ifdef MPI |
82 |
|
|
call mpi_barrier(mpi_comm_world,mpi_err) |
83 |
|
|
#endif |
84 |
|
|
end do dynamics |
85 |
|
|
|
86 |
|
|
end subroutine nose_hoover_dynamics |
87 |
|
|
|
88 |
|
|
|
89 |
|
|
subroutine hoovera() |
90 |
|
|
use velocity, ONLY : remove_cm_momentum,remove_angular_momentum |
91 |
|
|
|
92 |
|
|
|
93 |
|
|
! ** CONSTANT-TEMPERATURE MOLECULAR DYNAMICS USING CONSTRAINT. ** |
94 |
|
|
! ** ** |
95 |
|
|
! ** REFERENCES: ** |
96 |
|
|
! ** ** |
97 |
|
|
! ** HOOVER, LADD, AND MORAN, PHYS REV LETT 48, 1818, 1982. ** |
98 |
|
|
! ** EVANS, J CHEM PHYS 78, 3297, 1983. ** |
99 |
|
|
! ** BROWN AND CLARKE, MOL PHYS, 51, 1243, 1984. ** |
100 |
|
|
! ** ** |
101 |
|
|
! ** ROUTINES SUPPLIED: ** |
102 |
|
|
! ** ** |
103 |
|
|
! ** SUBROUTINE HOOVERA ( ) ** |
104 |
|
|
! ** FIRST PART OF MOVE WITH VELOCITY CONSTRAINTS APPLIED. ** |
105 |
|
|
! ** SUBROUTINE HOOVERB ( DT ) ** |
106 |
|
|
! ** SECOND PART OF MOVE. ** |
107 |
|
|
|
108 |
|
|
|
109 |
|
|
|
110 |
|
|
integer :: i, dim |
111 |
|
|
real( kind = DP) :: ke, v2, instatemp, chi |
112 |
|
|
real( kind = DP) :: kebar |
113 |
|
|
|
114 |
|
|
real( kind = DP ), dimension(ndim) :: v_dim_i |
115 |
|
|
|
116 |
|
|
real(kind=DP) :: ke_local |
117 |
|
|
|
118 |
|
|
! units for time are femtosec (10^-15 sec) |
119 |
|
|
! units for distance are angstroms (10^-10 m) |
120 |
|
|
! units for velocity are angstroms femtosec^-1 |
121 |
|
|
! units for mass are a.m.u. (1.661 * 10^-27 kg) |
122 |
|
|
! units for force are kcal mol^-1 angstrom^-1 |
123 |
|
|
! |
124 |
|
|
! converter will put the final terms into angstroms. |
125 |
|
|
! or angstrom/femtosecond. |
126 |
|
|
|
127 |
|
|
! converter = 1.0d0/2.390664d3 |
128 |
|
|
|
129 |
|
|
|
130 |
|
|
! first calculate the temperature using the unconstrained velocities: |
131 |
|
|
|
132 |
|
|
ke = 0.0E0_DP |
133 |
|
|
ke_local = 0.0E0_DP |
134 |
|
|
|
135 |
|
|
DO i = 1, nlocal |
136 |
|
|
do dim = 1, ndim |
137 |
|
|
v_dim_i(dim) = v(dim,i) + (dt2*f(dim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
138 |
|
|
ke_local = ke_local + mass(i)* & |
139 |
|
|
(v_dim_i(dim)*v_dim_i(dim)) |
140 |
|
|
enddo |
141 |
|
|
enddo |
142 |
|
|
#ifdef MPI |
143 |
|
|
call mpi_allreduce(ke_local,ke,1,mpi_double_precision, & |
144 |
|
|
mpi_sum,mpi_comm_world,mpi_err) |
145 |
|
|
#else |
146 |
|
|
ke = ke_local |
147 |
|
|
#endif |
148 |
|
|
|
149 |
|
|
|
150 |
|
|
ke = ke * AMUA2FS2_TO_KCALMOL |
151 |
|
|
kebar = kcal_to_kj*1000.0e0_DP*ke / real(3*(natoms - 1),DP) / NAVOGADRO |
152 |
|
|
instatemp = kebar/k_boltz |
153 |
|
|
|
154 |
|
|
!write(*,*) 'before, insta = ', instatemp |
155 |
|
|
|
156 |
|
|
chi = sqrt(target_temp / instatemp) |
157 |
|
|
|
158 |
|
|
! ** CALCULATE THE CONSTRAINED VELOCITIES AT TIME T+DT/2 ** |
159 |
|
|
|
160 |
|
|
do i = 1, nlocal |
161 |
|
|
do dim = 1, ndim |
162 |
|
|
v(dim,i) = v(dim,i) * (2.0E0_DP * chi - 1.0E0_DP) + & |
163 |
|
|
chi*dt*(f(dim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
164 |
|
|
enddo |
165 |
|
|
enddo |
166 |
|
|
|
167 |
|
|
if (sim_type == "cluster") then |
168 |
|
|
call remove_cm_momentum() |
169 |
|
|
call remove_angular_momentum() |
170 |
|
|
endif |
171 |
|
|
|
172 |
|
|
end subroutine hoovera |
173 |
|
|
|
174 |
|
|
|
175 |
|
|
|
176 |
|
|
subroutine hooverb( ) |
177 |
|
|
|
178 |
|
|
|
179 |
|
|
! ******************************************************************* |
180 |
|
|
! ** SECOND PART OF THE CONSTANT TEMPERATURE ALGORITHM ** |
181 |
|
|
! ** ** |
182 |
|
|
! ** THIS ADVANCES THE POSITIONS FROM T TO T + DT. ** |
183 |
|
|
! ******************************************************************* |
184 |
|
|
|
185 |
|
|
|
186 |
|
|
integer i, dim |
187 |
|
|
! units for time are femtosec (10^-15 sec) |
188 |
|
|
! units for distance are angstroms (10^-10 m) |
189 |
|
|
! units for velocity are angstroms femtosec^-1 |
190 |
|
|
! units for mass are a.m.u. (1.661 * 10^-27 kg) |
191 |
|
|
! units for force are kcal mol^-1 angstrom^-1 |
192 |
|
|
|
193 |
|
|
do i = 1, nlocal |
194 |
|
|
do dim = 1, ndim |
195 |
|
|
q(dim,i) = q(dim,i) + dt * v(dim,i) |
196 |
|
|
enddo |
197 |
|
|
enddo |
198 |
|
|
|
199 |
|
|
end subroutine hooverb |
200 |
|
|
|
201 |
|
|
|
202 |
|
|
end module dynamics_nose_hoover |