forked from lgcrego/Dynemol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
VDOS_m.f
369 lines (270 loc) · 9.34 KB
/
VDOS_m.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
!
! Calculates Vibrational spectral density and Density of states
! via the Fourier transform of the velocity autocorrelation function
!
! Alberto Torres
!
! input module: use and set in parameters_MM.f
module VDOS_input
! Sample length: Nr. of steps in each VACF sample: Nr. of v(0)'s
integer :: VDOS_Nsteps_per_sample = 10000 ! 0 -> turn off
end module VDOS_input
module VDOS_m
use constants_m , only : low_prec
use MD_read_m , only : MM, atom
use MM_input , only : species
use constants_m , only : twopi, h_bar, pico_2_sec
use parameters_m , only : n_t, t_i, t_f, restart, DRIVER
use VDOS_input , Lsample => VDOS_Nsteps_per_sample
public :: VDOS_init, VDOS_Correlation, VDOS_end
private
! module types ...
type velocity
real*8 :: vel(3)
end type velocity
! module parameters ...
integer, parameter :: in = 80, out = 81, out2 = 82
! module variables ...
integer :: Nres, & ! Number of residues (:= MM % N_of_species)
Natoms, &
Nsamples, &
isample, & ! sample index = [1,Nsamples] (this is the number of different 'v0')
frame_in_s ! frame index inside isample-th sample
type(velocity), allocatable :: v0(:), mw_v0(:)
real*8, allocatable :: v0_norm(:), v0_norm_mw(:)
real*8, allocatable :: VACF(:,:), VACF_mw(:,:) ! Velocity Auto Correlation Function
contains
! macro to expand dot product
#define VEC3_DOT_PROD(v,u) ( v(1)*u(1) + v(2)*u(2) + v(3)*u(3) )
! Parallelize if Natoms > SERIAL_THRESHOLD
#define SERIAL_THRESHOLD 500
!================
subroutine Get_v0
! Get initial velocity
!================
implicit none
! local variables
integer :: i, nr
v0_norm = 0.d0
v0_norm_mw = 0.d0
!$omp parallel do private(i,nr) default(shared) reduction(+:v0_norm,v0_norm_mw) if(Natoms > SERIAL_THRESHOLD)
do i = 1, Natoms
nr = atom(i) % nr
v0(i) % vel = atom(i) % vel + low_prec
v0_norm(nr) = v0_norm(nr) + VEC3_DOT_PROD( v0(i)%vel, v0(i)%vel )
mw_v0(i) % vel = ( v0(i) % vel )*( atom(i) % mass )
v0_norm_mw(nr) = v0_norm_mw(nr) + VEC3_DOT_PROD( mw_v0(i)%vel, v0(i)%vel )
end do
!$omp end parallel do
v0_norm(0) = sum( v0_norm (1:Nres) )
v0_norm_mw(0) = sum( v0_norm_mw(1:Nres) )
end subroutine Get_v0
!
!
!
!===================
subroutine VDOS_init
!===================
implicit none
! local variables
integer :: i
if (DRIVER == "q_dynamics" .or. &
DRIVER == "avrg_confgs" .or. &
DRIVER == "Genetic_Alg" .or. &
DRIVER == "diagnostic") Lsample = 0 ! turn off
! check Lsample's values
if (Lsample == 0) then
return
elseif (Lsample > n_t) then ! can't be greater than the number of steps
Lsample = n_t
elseif (mod(n_t,Lsample)/=0) then ! should be a perfect divisor of n_t: no fractional Nsamples
do while (mod(n_t,Lsample)/=0)
Lsample = Lsample + 1
end do
end if
Nres = MM % N_of_species
Natoms = MM % N_of_atoms
Nsamples = n_t / Lsample ! VACF will be averaged over "Nsamples" samples
write(*,*)
write(*,'(a)') " Velocity Autocorrelation: VDOS_init:"
write(*,'(a,i6)') " Total nr. of steps =", n_t
write(*,'(a,i6)') " Sample length (time steps) =", Lsample
write(*,'(a,i3)') " Nr. samples =", Nsamples
write(*,'(a,i3)') " Nr. residues =", Nres
write(*,*)
allocate( v0(Natoms), v0_norm (0:Nres) )
allocate( mw_v0(Natoms), v0_norm_mw(0:Nres) )
allocate( VACF( Lsample, 0:Nres ), source = 0.d0 )
allocate( VACF_mw( Lsample, 0:Nres ), source = 0.d0 )
if (restart) then
call VDOS_restart
return
end if
call Get_v0
isample = 1 ! sample index = [1,Nsamples] (this is the number of different 'v0's)
frame_in_s = 1 ! frame index inside isample-th sample
!------------------
! save v0 for restart (binary format)
open(out, file='VDOS.restart', status='unknown', form='unformatted')
write(out) v0_norm, v0_norm_mw, v0
close(out)
open(out, file='DOS.trunk/VACF.dat', status='unknown')
write(out,'(a)',advance='no') "# frame VACF: total "
do i = 1, Nres
write(out,'(a)',advance='no') species(i)%residue // ' '
end do
write(out,'(a)',advance='no') "MW-VACF: tot "
do i = 1, Nres
write(out,'(a)',advance='no') species(i)%residue // ' '
end do
write(out,*)
end subroutine VDOS_init
!
!
!
!======================
subroutine VDOS_restart
!======================
implicit none
open(in, file='VDOS.restart', status='old', form='unformatted', access='sequential', action='read')
read(in) v0_norm, v0_norm_mw, v0
close(in)
open(out, file='DOS.trunk/VACF.dat', status='old', access='append')
end subroutine VDOS_restart
!
!
!
!===================================
subroutine VDOS_Correlation( frame )
!===================================
implicit none
! inputs
integer, intent(in) :: frame ! global time step index
! local variables
integer :: i, j
real*8, allocatable :: auto_corr(:), auto_corr_mw(:) ! mw -> mass-weighted
if (Lsample == 0) return
allocate( auto_corr (0:Nres), source = 0.d0)
allocate( auto_corr_mw(0:Nres), source = 0.d0)
! Calculate <v(0).v(t)>
!$omp parallel do private(i,j) default(shared) reduction(+:auto_corr,auto_corr_mw) if(Natoms > SERIAL_THRESHOLD)
do i = 1, Natoms
j = atom(i) % nr
auto_corr(j) = auto_corr(j) + VEC3_DOT_PROD( v0(i)%vel, atom(i)%vel ) ! auto-correlalation
auto_corr_mw(j) = auto_corr_mw(j) + VEC3_DOT_PROD( mw_v0(i)%vel, atom(i)%vel ) ! mass-weighted auto-correlalation
end do
!$omp end parallel do
auto_corr (0) = sum( auto_corr (1:Nres) )
auto_corr_mw(0) = sum( auto_corr_mw(1:Nres) )
VACF (frame_in_s,:) = VACF (frame_in_s,:) + auto_corr (:) / v0_norm (:) ! accumulate normalized VACF
VACF_mw(frame_in_s,:) = VACF_mw(frame_in_s,:) + auto_corr_mw(:) / v0_norm_mw(:) ! accumulate normalized MW-VACF
write(out,'(i7)',advance='no') frame_in_s
do i = 0, Nres
write(out,'(es15.6)',advance='no') VACF(frame_in_s,i)/isample
end do
do i = 0, Nres
write(out,'(es15.6)',advance='no') VACF_mw(frame_in_s,i)/isample
end do
write(out,*)
if( mod(frame,Lsample) == 0 ) then
call Get_v0
isample = isample + 1
frame_in_s = 1
write(out,*)
return
end if
frame_in_s = frame_in_s + 1
end subroutine VDOS_Correlation
!
!
!
!==============================
subroutine VDOS_end
!==============================
implicit none
include 'fftw3.f'
! local variables
integer :: i, j
integer*8 :: plan
real*8 :: T, norm
real*8, parameter :: Hz_to_eV = twopi*pico_2_sec*h_bar
real*8, allocatable :: Fvacf(:,:) ! Fourier transform of the Auto Correlation
real*8, allocatable :: VDOS(:,:) ! Vibrational DOS
real*8, allocatable :: N_flex_atoms(:)
if (Lsample == 0) return
! release resources not needed anymore
deallocate(v0, mw_v0)
close(out) ! close VACF.dat
! count number of flexible atoms in each residue (for proprer normalization)
allocate( N_flex_atoms(0:Nres), source = 0.d0 )
do i = 1, Natoms
j = atom(i)%nr
if (atom(i)%flex) N_flex_atoms(j) = N_flex_atoms(j) + 1
end do
N_flex_atoms(0) = sum(N_flex_atoms(1:Nres)) ! total
!======= Vibrational Spectral Density ========
allocate( Fvacf(Lsample,0:Nres) )
VACF = VACF/Nsamples ! normalize average over samples
! using Real-Even symmetry: cosine fourier transform (refer to FFTW manual)
call dfftw_plan_r2r_1d( plan, Lsample, VACF(:,0), Fvacf(:,0), FFTW_REDFT00, FFTW_ESTIMATE )
do i = 0, Nres
call dfftw_execute_r2r( plan, VACF(:,i), Fvacf(:,i) )
end do
allocate( VDOS(Lsample,0:Nres) ) ! here, this is really the VSD (spectral density)
VDOS(:,:) = Fvacf(:,:)**2
! normalize
do i = 0, Nres
norm = (3*N_flex_atoms(i))/sum(VDOS(:,i))
VDOS(:,i) = VDOS(:,i)*norm
end do
deallocate( VACF ) ! not needed anymore
! Sample's Period (simulation time of each sample)
T = 2*(t_f - t_i)*pico_2_sec/Nsamples ! in seconds (factor 2 is to account for the real-even symmetry used above)
!------------------
! write to VSD-*.dat
do j = 0, Nres
if (j==0) then
open(out, file='DOS.trunk/VSD-total.dat', status='unknown')
else
open(out, file='DOS.trunk/VSD-'//trim(species(j)%residue)//'.dat', status='unknown')
end if
write(out,'(a)') "# E (eV) spectral density "
do i = 1, Lsample
write(out, '(2es)') ((i-1)/T)*Hz_to_eV, VDOS(i,j)
end do
close(out)
end do
!======== Vibrational Density of States ========
VACF_mw = VACF_mw/Nsamples ! normalize average over samples
do i = 0, Nres
call dfftw_execute_r2r( plan, VACF_mw(:,i), Fvacf(:,i) )
end do
call dfftw_destroy_plan( plan )
VDOS(:,:) = Fvacf(:,:)**2
! normalize
do i = 0, Nres
norm = (3*N_flex_atoms(i))/sum(VDOS(:,i))
VDOS(:,i) = VDOS(:,i)*norm
end do
!------------------
! write to VDOS-*.dat
do j = 0, Nres
if (j==0) then
open(out, file='DOS.trunk/VDOS-total.dat', status='unknown')
else
open(out, file='DOS.trunk/VDOS-'//trim(species(j)%residue)//'.dat', status='unknown')
end if
write(out,'(a)') "# E (eV) VDOS"
do i = 1, Lsample
write(out, '(2es)') ((i-1)/T)*Hz_to_eV, VDOS(i,j)
end do
close(out)
end do
! clean up
deallocate( VDOS, Fvacf, N_flex_atoms )
call system('rm -f VDOS.restart') ! if the calculation finishes, there's no more need of this file
end subroutine VDOS_end
!
!
!
end module VDOS_m