forked from lgcrego/Dynemol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
NVE.f
147 lines (116 loc) · 3.49 KB
/
NVE.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
module NVE_m
use constants_m
use syst ! using all syst
use MD_read_m , only: MM , atom , molecule
use VV_Parent , only: VV
public :: NVE
private
type, extends(VV) :: NVE
contains
procedure :: VV1
procedure :: VV2
end type NVE
interface NVE
module procedure constructor
end interface
contains
!
!
!
!===================================
function constructor() result( me )
!===================================
implicit none
type(NVE) :: 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(NVE) , intent(inout) :: me
real*8 , intent(in) :: dt
! local variables ...
real*8 :: ai(3)
real*8 :: massa , dt_HALF , dt2_HALF
integer :: i
dt_HALF = dt / two
dt2_HALF = dt_HALF * dt
! 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) % xyz = atom(i) % xyz + ( atom(i) % vel*dt + dt2_HALF*ai ) * mts_2_Angs
atom(i) % vel = atom(i) % vel + dt_HALF*ai
end if
end do
end subroutine VV1
!
!
!
!=========================
subroutine VV2( me , dt )
!=========================
implicit none
class(NVE) , 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
integer :: i , j , j1 , j2 , nresid
dt_half = dt / two
sumtemp = D_zero
! VV2 ...
select case (me % thermostat_type)
case (0:1) ! <== molecular ...
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
sumtemp = sumtemp + molecule(i) % mass * sum( V_CM * V_CM )
end do
case (2:) ! <== atomic ...
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
sumtemp = sumtemp + factor * sum( V_atomic * V_atomic )
end if
end do
end select
! 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
end subroutine VV2
!
!
!
end module NVE_m