forked from lgcrego/Dynemol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
NoseHoover.f
198 lines (162 loc) · 5.59 KB
/
NoseHoover.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
module Nose_Hoover_m
use constants_m
use syst ! using all syst
use MD_read_m , only: MM , atom , molecule , species
use VV_Parent , only: VV
use Berendsen_m , only: Barostat
public :: Nose_Hoover
private
type, extends(VV) :: Nose_Hoover
contains
procedure :: VV1
procedure :: VV2
end type Nose_Hoover
interface Nose_Hoover
module procedure constructor
end interface
! module variables ...
real*8 :: Csi, qmass, sigma
contains
!
!
!
!===================================
function constructor() result( me )
!===================================
implicit none
type(Nose_Hoover) :: me
!local variable ...
! select atomic or molecular kinetic energy to calculate the temperature ...
me % thermostat_type = NINT( float(maxval(molecule%N_of_atoms)) / float(MM%N_of_molecules) )
end function constructor
!
!
!
!========================
subroutine VV1( me , dt )
!========================
implicit none
class(Nose_Hoover) , intent(inout) :: me
real*8 , intent(in) :: dt
! local variables ...
real*8 :: ai(3) , tmp(3) , V_CM(3)
real*8 :: massa , dt_HALF , dt2_HALF
integer :: i , j , j1 , j2 , nresid
dt_HALF = dt / two
dt2_HALF = dt_HALF * dt
! calculation of the kinetic energy and thermostat-related things ...
me % kinetic = D_zero
select case (me % thermostat_type)
case (0:1) ! <== molecular ...
sigma = bath_T*boltz*THREE*real( count(molecule%flex) )*HALF
! molecular kinetic energy at time t ...
do i = 1 , MM % N_of_molecules
tmp = D_zero
nresid = molecule(i) % nr
j1 = sum(molecule(1:nresid-1) % N_of_atoms) + 1
j2 = sum(molecule(1:nresid) % N_of_atoms)
do j = j1 , j2
if ( atom(j) % flex ) then
massa = mol / atom(j) % mass
tmp = tmp + atom(j) % vel / massa
end if
end do
V_CM = tmp / molecule(i) % mass
me % kinetic = me % kinetic + molecule(i) % mass * sum( V_CM * V_CM ) * half
end do
case (2:) ! <== atomic ...
sigma = bath_T*boltz*THREE*real(count(atom%flex))*HALF
! atomic kinetic energy at time t ...
do i = 1 , MM % N_of_atoms
me % kinetic = me%kinetic + imol*atom(i)%mass*sum( atom(i) % vel(:) * atom(i) % vel(:) ) * half
end do
end select
qmass = TWO*sigma*(talt*pico_2_sec)**2
Csi = Csi + dt*( me%kinetic - sigma )/qmass
! VV1 ...
do i = 1 , MM % N_of_atoms
if( atom(i) % flex ) then
massa = mol / atom(i) % mass
ai = atom(i) % ftotal * massa
atom(i) % vel = atom(i) % vel * ( D_one - dt_HALF*Csi )
atom(i) % vel = atom(i) % vel + dt_HALF*ai
atom(i) % xyz = atom(i) % xyz + ( atom(i) % vel*dt ) * mts_2_Angs
end if
end do
end subroutine VV1
!
!
!
!=========================
subroutine VV2( me , dt )
!=========================
implicit none
class(Nose_Hoover) , intent(inout) :: me
real*8 , intent(in) :: dt
! local variables ...
real*8 :: tmp(3) , V_CM(3) , V_atomic(3)
real*8 :: massa , factor , sumtemp , dt_half , volume
integer :: i , j , j1 , j2 , nresid
dt_half = dt / two
sumtemp = D_zero
me % kinetic = D_zero
! VV2 ...
select case (me % thermostat_type)
case (0:1) ! <== molecular ...
! First, velocity and kinetic energy at (t + dt) ...
do i = 1 , MM % N_of_molecules
tmp = D_zero
nresid = molecule(i) % nr
j1 = sum(molecule(1:nresid-1) % N_of_atoms) + 1
j2 = sum(molecule(1:nresid) % N_of_atoms)
do j = j1 , j2
if ( atom(j) % flex ) then
massa = mol / atom(j) % mass
atom(j) % vel = atom(j) % vel + dt_half * atom(j) % ftotal * massa
tmp = tmp + atom(j) % vel / massa
end if
end do
V_CM = tmp / molecule(i) % mass
me % kinetic = me % kinetic + molecule(i) % mass * sum( V_CM * V_CM ) * half
sumtemp = sumtemp + molecule(i) % mass * sum( V_CM * V_CM )
end do
case (2:) ! <== atomic ...
! First, velocity and kinetic energy at (t + dt) ...
V_atomic = D_zero
do i = 1 , MM % N_of_atoms
if( atom(i) % flex ) then
massa = mol / atom(i) % mass
atom(i) % vel = atom(i) % vel + dt_half * atom(i) % ftotal * massa
V_atomic = atom(i) % vel
factor = imol * atom(i) % mass
me % kinetic = me % kinetic + factor * sum( V_atomic * V_atomic ) * half
sumtemp = sumtemp + factor * sum( V_atomic * V_atomic )
end if
end do
end select
Csi = Csi + dt*( me%kinetic - sigma )/qmass
do i = 1 , MM % N_of_atoms
atom(i) % vel = atom(i) % vel * ( D_one - dt_half*Csi )
end do
! instantaneous temperature of the system after contact with thermostat ...
select case (me % thermostat_type)
case (0:1) ! <== molecular ...
me % Temperature = sumtemp * iboltz / real( count(molecule%flex) )
case (2:) ! <atomic ...
me % Temperature = sumtemp * iboltz / real( count(atom%flex) )
end select
! calculation of the kinetic energy ...
me % Kinetic = D_zero
do i = 1 , MM % N_of_atoms
me % Kinetic = me % Kinetic + atom(i) % mass * sum( atom(i) % vel(:) * atom(i) % vel(:) ) * half
end do
me % Kinetic = me % Kinetic * micro / MM % N_of_Molecules
! calculation of pressure and density ...
massa = sum( species(:) % N_of_molecules * species(:) % mass )
volume = product( MM % box(:) * Angs_2_mts )
me % density = (massa / volume) * milli
end subroutine VV2
!
!
!
end module Nose_Hoover_m