-
Notifications
You must be signed in to change notification settings - Fork 0
/
18-molecule_npt_swap.f90
374 lines (282 loc) · 11 KB
/
18-molecule_npt_swap.f90
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
370
371
372
module mo_molecule_npt_swap
use mo_syst
use mo_config
use mo_list
use mo_molecule
use mo_extra_molecule
use mo_molecule_npt
use mo_math, only: swap
implicit none
!Assumption: 1) Integration method: prediction-correction; 2) Ensemble simulation: constraint method; 3) Soft particles; 4) Using neighbor list 5) The length unit is the smallest particle diameter
!procedure(abstract_extra_force), pointer :: tp_abstr_force => null() !the force type in this module shall be abstract_extra_force
type, extends(tpmolecule_npt), public :: tpmolecule_npt_swap
real(8) :: swap_prob
contains
!override
procedure, pass :: initial => initial_npt_swap
procedure, pass :: body => body_npt_swap
procedure, pass :: rescale_p => rescale_p_npt_swap
!extend
procedure, pass :: sweep_swap => sweep_swap
procedure, pass :: pair_swap => pair_swap
!inherit :: variables, clean_npt, predic_npt, correc_npt, constri_npt, rescale_t_npt
end type tpmolecule_npt_swap
contains
subroutine initial_npt_swap( this, tcon, tnb, tp_abstr_force, tset_molecule )
implicit none
! para list
class(tpmolecule_npt_swap), intent(inout) :: this
type(con_t), intent(inout) :: tcon
type(list_t), intent(inout) :: tnb
class(tpset_molecule), intent(in) :: tset_molecule
logical, external :: tp_abstr_force
!local
integer :: i
logical :: get_force
select type(tset_molecule)
type is (tpset_npt_swap)
associate( &
unit_v => this%unit_v, &
temper => this%temper, &
pre => this%pre, &
dt => this%dt, &
c1 => this%c1, &
c2 => this%c2, &
c3 => this%c3, &
coeff0 => this%coeff0, &
coeff2 => this%coeff2, &
coeff3 => this%coeff3, &
la_list => this%la_list,&
la1 => this%la1, &
la2 => this%la2, &
la3 => this%la3, &
lagr_t => this%lagr_t, &
lagr_p => this%lagr_p, &
natom => tcon%natom, &
radius => tcon%r, &
ra => tcon%ra, &
va => tcon%va, &
fa => tcon%fa, &
la => tcon%la &
)
if ( .not. allocated(this%ra1) ) allocate( this%ra1(free,natom) )
if ( .not. allocated(this%ra2) ) allocate( this%ra2(free,natom) )
if ( .not. allocated(this%ra3) ) allocate( this%ra3(free,natom) )
if ( .not. allocated(this%va1) ) allocate( this%va1(free,natom) )
if ( .not. allocated(this%va2) ) allocate( this%va2(free,natom) )
if ( .not. allocated(this%va3) ) allocate( this%va3(free,natom) )
this%swap_prob = tset_molecule%swap_prob
this%resc_prob = tset_molecule%resc_prob
temper = tset_molecule%temper
pre = tset_molecule%pre
dt = tset_molecule%dt
c1 = dt
c2 = c1 * dt / 2.d0
c3 = c2 * dt / 3.d0
coeff0 = this%gear0 * c1
coeff2 = this%gear2 * c1 / c2
coeff3 = this%gear3 * c1 / c3
!unit_v = sum( (2.d0*radius)**dble(free) ) / dble(natom)
unit_v = 1.d0
call full_make_list( tnb, tcon)
la_list = la
this%nlcut0 = tnb%nlcut
call random_number(va)
call rescale_t_npt( this, tcon )
call rescale_p_npt_swap( this, tcon, tnb, tp_abstr_force )
get_force = tp_abstr_force( tcon, tnb, 0 )
call constri_npt( this, tcon )
associate( &
ra1 => this%ra1, &
ra2 => this%ra2, &
ra3 => this%ra3, &
va1 => this%va1, &
va2 => this%va2, &
va3 => this%va3 &
)
ra1 = va + lagr_p * ra
ra2 = 0.d0
ra3 = 0.d0
va1 = fa - lagr_t * va
va2 = 0.d0
va3 = 0.d0
la1 = lagr_p * la
la2 = 0.d0
la3 = 0.d0
end associate
end associate
end select
end subroutine
subroutine body_npt_swap( this, tcon, tnb, tp_abstr_force, tout_molecule )
implicit none
! para list
class(tpmolecule_npt_swap), intent(inout) :: this
type(con_t), intent(inout) :: tcon
type(list_t), intent(inout) :: tnb
logical, external :: tp_abstr_force
type(tpout_molecule), intent(inout), optional :: tout_molecule
! local
integer :: i
logical :: get_force
real(8) :: prob, temp
call predic_npt( this, tcon )
call sweep_swap( this, tcon, tnb, tp_abstr_force )
temp = dble(free) * (tcon%la(1) - this%la_list(1))
if(temp < 0) then
tnb%nlcut = this%nlcut0 + temp
endif
if ( tnb%nlcut < 1.d-2 .or. check_list( tnb, tcon )) then
call full_make_list( tnb, tcon )
this%la_list = tcon%la
end if
get_force = tp_abstr_force( tcon, tnb, 0 )
call constri_npt( this, tcon )
call correc_npt( this, tcon )
if ( present( tout_molecule ) ) then
associate( &
unit_v => this%unit_v, &
natom => tcon%natom, &
radius => tcon%r, &
la => tcon%la, &
vili => ex_var%vili &
)
tout_molecule%phi = sqrt( pi**free ) / gamma( dble(free)/2.d0+1.d0 ) * sum( radius**free ) / product( la )
tout_molecule%T = 2.d0 * tcon%Ek / dble(free)
tout_molecule%press = ( tout_molecule%T * dble(natom) + vili ) * unit_v / product( la )
tout_molecule%Ek = tcon%Ek
tout_molecule%Ea = tcon%Ea
tout_molecule%Ev = tcon%Ek + tcon%Ea
end associate
end if
call random_number( prob )
if( prob < this%resc_prob ) then
call rescale_t_npt( this, tcon )
call rescale_p_npt_swap( this, tcon, tnb, tp_abstr_force )
end if
end subroutine
subroutine sweep_swap( this, tcon, tnb, tp_abstr_force )
implicit none
! para list
class(tpmolecule_npt_swap), intent(inout) :: this
type(con_t), intent(inout) :: tcon
type(list_t), intent(inout) :: tnb
logical, external :: tp_abstr_force
!local
integer :: n, ti, tj, tn
real(8) :: var, prob
associate( &
natom => tcon%natom, &
radius => tcon%r, &
swap_prob => this%swap_prob &
)
do n = 1, natom
call random_number(prob)
if(prob < swap_prob) then
call random_number(var)
ti = floor( var * dble(natom) ) + 1
call random_number(var)
tj = floor( var * dble(natom) ) + 1
if( ti .eq. tj ) cycle
if( pair_swap( this, tcon, tnb, tp_abstr_force, ti, tj ) ) then
if( radius(ti) > radius(tj) ) then
tn = ti
else
tn = tj
end if
call one_make_list( tnb, tcon, tn )
end if
end if
enddo
end associate
end subroutine
logical function pair_swap( this, tcon, tnb, tp_abstr_force, ti, tj )
implicit none
! para list
class(tpmolecule_npt_swap), intent(inout) :: this
type(con_t), intent(inout) :: tcon
type(list_t), intent(inout) :: tnb
logical, external :: tp_abstr_force
integer, intent(in) :: ti, tj
!local
real(8) :: prob, temp
real(8) :: ea0, ea1
logical :: get_potential
associate( &
temper => this%temper, &
radius => tcon%r, &
dea => ex_var%dea &
)
pair_swap = .false.
temp = abs( radius(ti) - radius(tj) )
if(temp > tnb%nlcut) then
go to 100
end if
ea0 = 0.d0
get_potential = tp_abstr_force( tcon, tnb, ti )
ea0 = ea0 + dea
get_potential = tp_abstr_force( tcon, tnb, tj )
ea0 = ea0 + dea
call swap( radius(ti), radius(tj) )
ea1 = 0.d0
get_potential = tp_abstr_force( tcon, tnb, ti )
ea1 = ea1 + dea
get_potential = tp_abstr_force( tcon, tnb, tj )
ea1 = ea1 + dea
if( ea1 < ea0 ) then
pair_swap = .true.
go to 100
end if
prob = exp( (ea0 - ea1)/temper )
call random_number(temp)
if( temp < prob ) then
pair_swap = .true.
go to 100
else
call swap( radius(ti), radius(tj) )
end if
100 end associate
return
end function
subroutine rescale_p_npt_swap( this, tcon, tnb, tp_abstr_force )
implicit none
! para list
class(tpmolecule_npt_swap), intent(inout) :: this
type(con_t), intent(inout) :: tcon
type(list_t), intent(inout) :: tnb
logical, external :: tp_abstr_force
!local
real(8) :: pvar, ttemper, temp, factor
logical :: get_force
associate( &
unit_v => this%unit_v, &
pre => this%pre, &
la_list => this%la_list, &
natom => tcon%natom, &
ra => tcon%ra, &
la => tcon%la, &
chixi => ex_var%chixi, &
vili => ex_var%vili &
)
ttemper = sum( tcon%va**2.d0 ) / dble(free)
call full_make_list( tnb, tcon )
la_list = la
get_force = tp_abstr_force( tcon, tnb, 0 )
pvar = ( ttemper + vili ) * unit_v / product( la )
do while( abs( (pre-pvar)/pre ) > 1.d-2)
factor = 1.d0 + (pvar - pre) / (pvar + chixi/dble(free**2)/product(la) ) / dble(free)
la = la * factor
ra = ra * factor
temp = dble(free) * (tcon%la(1) - this%la_list(1))
if(temp < 0) then
tnb%nlcut = this%nlcut0 + temp
endif
if ( tnb%nlcut < 1.d-2 ) then
call full_make_list( tnb, tcon )
this%la_list = tcon%la
end if
get_force = tp_abstr_force( tcon, tnb, 0 )
pvar = ( ttemper + vili ) * unit_v / product( la )
end do
end associate
end subroutine
end module