1 |
chuckv |
4 |
!=========================================================================== |
2 |
|
|
! |
3 |
|
|
!=========================================================================== |
4 |
|
|
! . Subroutines: |
5 |
|
|
! LANGEVIN_VERLET_DYNAMICS molecular dynamics using langevin |
6 |
|
|
! dynamics |
7 |
|
|
!=========================================================================== |
8 |
|
|
|
9 |
|
|
|
10 |
|
|
module dynamics_langevin_verlet |
11 |
|
|
use definitions, ONLY : DP,NDIM |
12 |
|
|
use constants, ONLY : AMU_TO_KG,K_BOLTZ,MS_TO_AFS, & |
13 |
|
|
KCALMOL_TO_AMUA2FS2, PI |
14 |
|
|
use simulation |
15 |
|
|
use parameter |
16 |
|
|
use force_module, ONLY : calc_forces |
17 |
|
|
use second_deriv |
18 |
|
|
use force_utilities, ONLY : check |
19 |
|
|
use status |
20 |
|
|
! use velocity, ONLY : scale_velocities |
21 |
|
|
use dynamics_utilities, ONLY : evolve_time, init_dynamics_loop, & |
22 |
|
|
sim_status,DT2,DTSQRT,DTSQ2, init_status |
23 |
|
|
use velocity, ONLY : calc_temp |
24 |
|
|
#ifdef MPI |
25 |
|
|
use mpi_module, ONLY : mpi_comm_world,mpi_err |
26 |
|
|
#endif |
27 |
|
|
|
28 |
|
|
implicit none |
29 |
|
|
PRIVATE |
30 |
|
|
SAVE |
31 |
|
|
|
32 |
|
|
! . Constants common to langevin dynamcics |
33 |
|
|
|
34 |
|
|
real ( kind = dp ) :: C0,C1,C2 |
35 |
|
|
real ( kind = dp ), allocatable, dimension(:) :: CRV1,CRV2 |
36 |
|
|
real ( kind = dp ), allocatable, dimension(:) :: FACR1,FACR2 |
37 |
|
|
real ( kind = dp ), allocatable, dimension(:) :: FACV1,FACV2,FACV3 |
38 |
|
|
|
39 |
|
|
real ( kind = dp ), allocatable, dimension(:) :: SDR,SDV |
40 |
|
|
|
41 |
|
|
|
42 |
|
|
! . langevin arrays |
43 |
|
|
real ( kind = dp ), allocatable, dimension(:,:) :: v_lang |
44 |
|
|
real ( kind = dp ), allocatable, dimension(:) :: MASFAC |
45 |
|
|
|
46 |
|
|
type (sim_status) :: langevin_status |
47 |
|
|
|
48 |
|
|
public :: langevin_verlet_dynamics |
49 |
|
|
|
50 |
|
|
|
51 |
|
|
|
52 |
|
|
|
53 |
|
|
contains |
54 |
|
|
|
55 |
|
|
! . subroutine init_langevin sets up common |
56 |
|
|
! constants for dynamics |
57 |
|
|
subroutine init_langevin |
58 |
|
|
use model_module, ONLY : get_max_atype |
59 |
|
|
integer :: n_atypes_known |
60 |
|
|
character(len=80) :: msg |
61 |
|
|
real( kind = dp ) :: fact |
62 |
|
|
real( kind = dp ) :: this_eta |
63 |
|
|
! . allocate nlocal arrays |
64 |
|
|
allocate(v_lang(ndim,nlocal)) |
65 |
|
|
allocate(masfac(nlocal)) |
66 |
|
|
! . find out how many atypes the simulation |
67 |
|
|
! . knows about |
68 |
|
|
|
69 |
|
|
n_atypes_known = get_max_atype() |
70 |
|
|
allocate(SDR(n_atypes_known)) |
71 |
|
|
allocate(SDV(n_atypes_known)) |
72 |
|
|
allocate(CRV1(n_atypes_known)) |
73 |
|
|
allocate(CRV2(n_atypes_known)) |
74 |
|
|
allocate(FACR1(n_atypes_known)) |
75 |
|
|
allocate(FACR2(n_atypes_known)) |
76 |
|
|
allocate(FACV1(n_atypes_known)) |
77 |
|
|
allocate(FACV2(n_atypes_known)) |
78 |
|
|
allocate(FACV3(n_atypes_known)) |
79 |
|
|
|
80 |
|
|
|
81 |
|
|
this_eta = eta |
82 |
|
|
call calc_gamma_factors(this_eta) |
83 |
|
|
|
84 |
|
|
! . Temperature and mass factors. |
85 |
|
|
FACT = ( K_BOLTZ * bath_temp / AMU_TO_KG ) |
86 |
|
|
!!!!!DANGER DANGER WILL |
87 |
|
|
! WHERE ( .NOT. ATMFIX ) MASFAC = MS_TO_AFS * SQRT ( FACT / mass ) ! A fs^-1. |
88 |
|
|
MASFAC = MS_TO_AFS * SQRT ( FACT / mass ) ! A fs^-1. |
89 |
|
|
|
90 |
|
|
write(msg,*) "Visocity of solvent is: ", eta |
91 |
|
|
call info("INIT_LANGEVIN",msg) |
92 |
|
|
|
93 |
|
|
|
94 |
|
|
end subroutine init_langevin |
95 |
|
|
|
96 |
|
|
|
97 |
|
|
subroutine langevin_verlet_dynamics() |
98 |
|
|
use langevin_cluster_list, ONLY : init_langevin_list,create_langevin_list |
99 |
|
|
|
100 |
|
|
|
101 |
|
|
logical :: update_nlist |
102 |
|
|
logical :: do_pot |
103 |
|
|
logical :: not_done |
104 |
|
|
logical :: nmflag |
105 |
|
|
|
106 |
|
|
integer :: step, i |
107 |
|
|
real(kind = dp) :: temp |
108 |
|
|
call init_status(langevin_status) |
109 |
|
|
call init_dynamics_loop(langevin_status) |
110 |
|
|
call init_langevin_list(nlocal) |
111 |
|
|
! Start the main simulation loop. |
112 |
|
|
|
113 |
|
|
do_pot = .true. |
114 |
|
|
not_done = .true. |
115 |
|
|
nmflag = .false. |
116 |
|
|
step = 0 |
117 |
|
|
!. initialize constants for dynamics |
118 |
|
|
call init_langevin() |
119 |
|
|
call create_langevin_list() |
120 |
|
|
dynamics: do |
121 |
|
|
if (.not.not_done) exit dynamics |
122 |
|
|
|
123 |
|
|
step = step + 1 |
124 |
|
|
|
125 |
|
|
|
126 |
|
|
! if (checktemp.and.(mod(step,check_temp_steps).eq.0)) then |
127 |
|
|
! call calc_temp(langevin_status%temp) |
128 |
|
|
! if (dabs(langevin_status%temp-target_temp).gt.therm_variance) then |
129 |
|
|
! call scale_velocities(target_temp) |
130 |
|
|
! end if |
131 |
|
|
! end if |
132 |
|
|
|
133 |
|
|
|
134 |
|
|
|
135 |
|
|
call check(update_nlist) |
136 |
|
|
if (update_nlist) then |
137 |
|
|
call create_langevin_list() |
138 |
|
|
endif |
139 |
|
|
|
140 |
|
|
call langevina() |
141 |
|
|
|
142 |
|
|
! . Add in random forces |
143 |
|
|
call RANDOM_FORCES |
144 |
|
|
|
145 |
|
|
if (do_pot) then |
146 |
|
|
call calc_forces(update_nlist,nmflag,pe = langevin_status%pot_e) |
147 |
|
|
else |
148 |
|
|
call calc_forces(update_nlist,nmflag) |
149 |
|
|
endif |
150 |
|
|
|
151 |
|
|
|
152 |
|
|
call langevinb() |
153 |
|
|
|
154 |
|
|
call evolve_time(step,langevin_status,do_pot,not_done) |
155 |
|
|
|
156 |
|
|
|
157 |
|
|
#ifdef MPI |
158 |
|
|
call mpi_barrier(mpi_comm_world,mpi_err) |
159 |
|
|
#endif |
160 |
|
|
end do dynamics |
161 |
|
|
|
162 |
|
|
deallocate(v_lang) |
163 |
|
|
deallocate(masfac) |
164 |
|
|
end subroutine langevin_verlet_dynamics |
165 |
|
|
|
166 |
|
|
|
167 |
|
|
! . First part of langevin dynamics |
168 |
|
|
! Advances positions from T to T + DT |
169 |
|
|
! Advances velocities from T to T + DT/2 |
170 |
|
|
!=================================== |
171 |
|
|
subroutine langevina() |
172 |
|
|
!=================================== |
173 |
|
|
use langevin_cluster_list, ONLY : langevin_list |
174 |
|
|
integer :: i |
175 |
|
|
! units for time are femtosec (10^-15 sec) |
176 |
|
|
! units for distance are angstroms (10^-10 m) |
177 |
|
|
! units for velocity are angstroms femtosec^-1 |
178 |
|
|
! units for mass are a.m.u. (1.661 * 10^-27 kg) |
179 |
|
|
! units for force are kcal mol^-1 angstrom^-1 |
180 |
|
|
! |
181 |
|
|
! converter will put the final terms into angstroms. |
182 |
|
|
! or angstrom/femtosecond. |
183 |
|
|
|
184 |
|
|
! Calculate new positions and velocities |
185 |
|
|
v_lang = 0.0_DP |
186 |
|
|
|
187 |
|
|
do i = 1, nlocal |
188 |
|
|
if (langevin_list(i)) then |
189 |
|
|
|
190 |
|
|
q(1:ndim,i) = q(1:ndim,i) + & |
191 |
|
|
FACR1(ident(i)) * v(1:ndim,i) + & |
192 |
|
|
(FACR2(ident(i)) * f(1:ndim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
193 |
|
|
v_lang(1:ndim,i) = FACV1(ident(i)) * v(1:ndim,i) + & |
194 |
|
|
(FACV2(ident(i)) * f(1:ndim,i)/mass(i)) * KCALMOL_TO_AMUA2FS2 |
195 |
|
|
|
196 |
|
|
else |
197 |
|
|
q(1:ndim,i) = q(1:ndim,i) + dt*v(1:ndim,i) & |
198 |
|
|
+ (dtsq2*f(1:ndim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
199 |
|
|
v(1:ndim,i) = v(1:ndim,i) & |
200 |
|
|
+ (dt2*f(1:ndim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
201 |
|
|
|
202 |
|
|
endif |
203 |
|
|
|
204 |
|
|
|
205 |
|
|
end do |
206 |
|
|
|
207 |
|
|
end subroutine langevina |
208 |
|
|
|
209 |
|
|
|
210 |
|
|
!======================== |
211 |
|
|
subroutine langevinb( ) |
212 |
|
|
!========================= |
213 |
|
|
use langevin_cluster_list, ONLY : langevin_list |
214 |
|
|
integer :: i |
215 |
|
|
! units for time are femtosec (10^-15 sec) |
216 |
|
|
! units for distance are angstroms (10^-10 m) |
217 |
|
|
! units for velocity are angstroms femtosec^-1 |
218 |
|
|
! units for mass are a.m.u. (1.661 * 10^-27 kg) |
219 |
|
|
! units for force are kcal mol^-1 angstrom^-1 |
220 |
|
|
|
221 |
|
|
do i = 1, nlocal |
222 |
|
|
if (langevin_list(i)) then |
223 |
|
|
v(1:ndim,i) = v_lang(1:ndim,i) & |
224 |
|
|
+ (FACV3(ident(i)) * f(1:ndim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
225 |
|
|
else ! standard movea |
226 |
|
|
v(1:ndim,i) = v(1:ndim,i) & |
227 |
|
|
+ (dt2*f(1:ndim,i)/mass(i))*KCALMOL_TO_AMUA2FS2 |
228 |
|
|
endif |
229 |
|
|
enddo |
230 |
|
|
|
231 |
|
|
end subroutine langevinb |
232 |
|
|
|
233 |
|
|
|
234 |
|
|
|
235 |
|
|
!----------------------- |
236 |
|
|
SUBROUTINE RANDOM_FORCES |
237 |
|
|
!----------------------- |
238 |
|
|
use langevin_cluster_list, ONLY : langevin_list |
239 |
|
|
use sprng_mod, ONLY : random_gauss |
240 |
|
|
! . Local scalars. |
241 |
|
|
INTEGER :: I, IATOM |
242 |
|
|
REAL ( KIND = DP ) :: RANR, RANV, R1, R2 |
243 |
|
|
|
244 |
|
|
! . Loop over the atoms. |
245 |
|
|
DO IATOM = 1,nlocal |
246 |
|
|
|
247 |
|
|
! . The atom is free. |
248 |
|
|
IF ( langevin_list(iatom) ) THEN |
249 |
|
|
|
250 |
|
|
! . Loop over the Cartesian components for the atom. |
251 |
|
|
DO I = 1,ndim |
252 |
|
|
|
253 |
|
|
! . Generate two Gaussian random numbers. |
254 |
|
|
R1 = RANDOM_GAUSS ( ) |
255 |
|
|
R2 = RANDOM_GAUSS ( ) |
256 |
|
|
|
257 |
|
|
! . Calculate the position and velocity changes. |
258 |
|
|
RANR = MASFAC(IATOM) * SDR(ident(iatom)) * R1 |
259 |
|
|
RANV = MASFAC(IATOM) * SDV(ident(iatom)) * & |
260 |
|
|
( CRV1(ident(iatom)) * R1 + CRV2(ident(iatom)) * R2 ) |
261 |
|
|
|
262 |
|
|
! . Add in the components. |
263 |
|
|
q(I,IATOM) = q(I,IATOM) + RANR |
264 |
|
|
v_lang(I,IATOM) = v_lang(I,IATOM) + RANV |
265 |
|
|
|
266 |
|
|
END DO |
267 |
|
|
END IF |
268 |
|
|
END DO |
269 |
|
|
|
270 |
|
|
END SUBROUTINE RANDOM_FORCES |
271 |
|
|
|
272 |
|
|
SUBROUTINE calc_gamma_factors(this_eta) |
273 |
|
|
use model_module, ONLY : atype_mass,atype_identifier |
274 |
|
|
integer :: atype, i |
275 |
|
|
character(len = 80) :: msg |
276 |
|
|
real( kind = dp ), intent(in) :: this_eta |
277 |
|
|
real( kind = dp ) :: gamma |
278 |
|
|
real( kind = dp ) :: FACT |
279 |
|
|
real( kind = dp ) :: C0 |
280 |
|
|
real( kind = dp ) :: C1 |
281 |
|
|
real( kind = dp ) :: C2 |
282 |
|
|
|
283 |
|
|
|
284 |
|
|
do i = 1, n_atypes_present |
285 |
|
|
gamma = 0.0_DP |
286 |
|
|
fact = 0.0_DP |
287 |
|
|
C0 = 0.0_DP |
288 |
|
|
C1 = 0.0_DP |
289 |
|
|
C2 = 0.0_DP |
290 |
|
|
|
291 |
|
|
atype = present_atypes(i) |
292 |
|
|
!. calculate Gamma in fs^-1 |
293 |
|
|
!. gamma = (6 * pi * atom sigma(angstroms) * viscosity(poise) * 100cm/m) * 1fs/10^-15sec / |
294 |
|
|
!. (atom mass(amu) * amu_to_kg * 1000g/kg) |
295 |
|
|
|
296 |
|
|
gamma = (6 * pi * get_langevin_sigma(atype) * this_eta ) / & |
297 |
|
|
(atype_mass(atype)*amu_to_kg*10) * 1E-15_DP |
298 |
|
|
write(msg,*) "Gamma for ",atype_identifier(atype), "is: ", gamma |
299 |
|
|
call info("CALC_GAMMA_FACTORS",msg) |
300 |
|
|
|
301 |
|
|
FACT = dt * gamma |
302 |
|
|
|
303 |
|
|
C0 = EXP ( - FACT ) |
304 |
|
|
C1 = ( 1.0_DP - C0 ) / FACT |
305 |
|
|
C2 = ( 1.0_DP - C1 ) / FACT |
306 |
|
|
|
307 |
|
|
! . Random force constants. |
308 |
|
|
SDR(atype) = SQRT ( dt * dt * & |
309 |
|
|
( 2.0_DP - ( 3.0_DP - 4.0_DP * C0 + C0 * C0 ) / FACT ) / FACT ) ! fs. |
310 |
|
|
SDV(atype) = SQRT ( ( 1.0_DP - C0 * C0 ) ) ! Dimensionless. |
311 |
|
|
CRV1(atype) = dt * ( 1.0_DP - C0 )**2 / ( FACT * SDR(atype) * SDV(atype) ) ! Dimensionless. |
312 |
|
|
CRV2(atype) = SQRT ( 1.0_DP - CRV1(atype) * CRV1(atype) ) ! Dimensionless. |
313 |
|
|
|
314 |
|
|
|
315 |
|
|
! . Calculate some integration constants. |
316 |
|
|
FACR1(atype) = C1 * dt ! fs. |
317 |
|
|
FACR2(atype) = C2 * dt * dt ! fs^2. |
318 |
|
|
FACV1(atype) = C0 ! Dimensionless. |
319 |
|
|
FACV2(atype) = ( C1 - C2 ) * dt ! fs. |
320 |
|
|
FACV3(atype) = C2 * dt ! fs. |
321 |
|
|
|
322 |
|
|
end do |
323 |
|
|
END SUBROUTINE calc_gamma_factors |
324 |
|
|
|
325 |
|
|
|
326 |
|
|
function get_langevin_sigma(this_atype) result(this_sigma) |
327 |
|
|
integer, intent(in) :: this_atype |
328 |
|
|
real( kind = dp ) :: this_sigma |
329 |
|
|
include 'headers/atom.h' |
330 |
|
|
|
331 |
|
|
select case (this_atype) |
332 |
|
|
|
333 |
|
|
case(Au_atom) |
334 |
|
|
this_sigma = 1.647E-10_DP |
335 |
|
|
case(Ag_atom) |
336 |
|
|
this_sigma = 1.574E-10_DP |
337 |
|
|
case(Pt_atom) |
338 |
|
|
this_sigma = 1.377E-10_DP |
339 |
|
|
case(Pb_atom) |
340 |
|
|
this_sigma = 2.148E-10_DP |
341 |
|
|
case(Cu_atom) |
342 |
|
|
this_sigma = 1.748E-10_DP |
343 |
|
|
case(Ni_atom) |
344 |
|
|
this_sigma = 1.417E-10_DP |
345 |
|
|
case(Pd_atom) |
346 |
|
|
this_sigma = 1.450E-10_DP |
347 |
|
|
case default |
348 |
|
|
call error("get_langevin_sigma","Unknown atom type") |
349 |
|
|
end select |
350 |
|
|
end function get_langevin_sigma |
351 |
|
|
end module dynamics_langevin_verlet |