forked from libAtoms/pymatnest
-
Notifications
You must be signed in to change notification settings - Fork 0
/
example_bead_spring_polymer_model.F90
373 lines (328 loc) · 13 KB
/
example_bead_spring_polymer_model.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
373
subroutine ll_init_model(N_params, params)
use example_bead_spring_polymer_params_mod
implicit none
integer :: N_params
double precision :: params(N_params)
! PARAMS: chain_length K_bond r0_bond K_angle epsilon sigma cutoff
if (N_params /= 7) then
print *, "example_bead_spring_polymer_model.F90 got 7 /= N_params ", N_params
print *, "Need: chain_length K_bond r0_bond K_angle epsilon sigma cutoff"
call exit(1)
endif
chain_length = int(params(1))
K_bond = params(2)
r0_bond = params(3)
K_angle = params(4)
epsilon = params(5)
sigma = params(6)
cutoff = params(7)
cutoff_sq = cutoff*cutoff
E_offset = (sigma/cutoff)**12 - (sigma/cutoff)**6
end subroutine ll_init_model
subroutine ll_init_config(N, Z, pos, cell, Emax)
implicit none
integer :: N
integer :: Z(N)
double precision :: pos(3,N), cell(3,3)
double precision :: Emax
return
end subroutine ll_init_config
double precision function ll_eval_energy(N, Z, pos, n_extra_data, extra_data, cell)
use example_mat_mod
use example_bead_spring_polymer_params_mod
implicit none
integer :: N
integer :: Z(N)
double precision :: pos(3,N), cell(3,3)
integer :: n_extra_data
double precision :: extra_data(n_extra_data, N)
integer :: i, j
double precision :: dr(3), dr_mag, dr_mag_sq, dr_l(3), dr_l0(3), pos_l(3,N), drr(3)
double precision :: cell_inv(3,3), E_term
integer :: dj1, dj2, dj3
integer :: n_images
double precision cell_height(3), v_norm_hat(3)
integer :: i_chain, j_chain, n_chains, i_bond
call matrix3x3_inverse(cell, cell_inv)
! into lattice coodinates
pos_l = matmul(cell_inv, pos)
if (n_extra_data == 1) extra_data = 0.0
do i=1, 3
v_norm_hat = cell(:,mod(i,3)+1) .cross. cell(:,mod(i+1,3)+1)
v_norm_hat = v_norm_hat / sqrt(sum(v_norm_hat**2))
cell_height(i) = abs(sum(v_norm_hat*cell(:,i)))
end do
n_images = ceiling(cutoff/minval(cell_height))
ll_eval_energy = 0.0
do i=1, N
i_chain = int((i-1)/chain_length)
do j=i, N
j_chain = int((j-1)/chain_length)
dr_l0 = pos_l(:,i)-pos_l(:,j)
dr_l0 = dr_l0 - floor(dr_l0+0.5)
do dj1=-n_images,n_images
dr_l(1) = dr_l0(1) + real(dj1, 8)
do dj2=-n_images,n_images
dr_l(2) = dr_l0(2) + real(dj2, 8)
do dj3=-n_images,n_images
dr_l(3) = dr_l0(3) + real(dj3, 8)
if (i == j .and. dj1 == 0 .and. dj2 == 0 .and. dj3 == 0) cycle ! no self-interaction
if (i_chain == j_chain .and. abs(i-j) == 1 .and. dj1 == 0 .and. dj2 == 0 .and. dj3 == 0) cycle ! bond excluded
dr(1) = cell(1,1) * dr_l(1) + cell(1,2) * dr_l(2) + cell(1,3) * dr_l(3) ! sum(cell(1,:)*dr_l)
dr(2) = cell(2,1) * dr_l(1) + cell(2,2) * dr_l(2) + cell(2,3) * dr_l(3) ! sum(cell(2,:)*dr_l)
dr(3) = cell(3,1) * dr_l(1) + cell(3,2) * dr_l(2) + cell(3,3) * dr_l(3) ! sum(cell(3,:)*dr_l)
dr_mag_sq = dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3) ! sum(dr*dr)
if (dr_mag_sq < cutoff_sq) then
dr_mag = sqrt(dr_mag_sq)
E_term = epsilon*(((sigma/dr_mag)**12 - (sigma/dr_mag)**6) - E_offset)
if (i == j) E_term = E_term * 0.5
ll_eval_energy = ll_eval_energy + E_term
if (n_extra_data == 1) then
extra_data(1,i) = extra_data(1,i) + 0.5*E_term
extra_data(1,j) = extra_data(1,j) + 0.5*E_term
endif
endif
end do
end do
end do
end do
end do
n_chains = N/chain_length
do i_chain=1, n_chains
! bonds
do i_bond=(i_chain-1)*chain_length+1, i_chain*chain_length-1
dr = pos(:,i_bond+1) - pos(:, i_bond)
dr_mag = sqrt(dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3))
ll_eval_energy = ll_eval_energy + K_bond * (dr_mag-r0_bond)*(dr_mag-r0_bond)
end do
! angles
do i_bond=(i_chain-1)*chain_length+2, i_chain*chain_length-1
dr = pos(:,i_bond+1) - pos(:, i_bond)
dr = dr / sqrt(dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3))
drr = pos(:,i_bond-1) - pos(:, i_bond)
drr = drr / sqrt(drr(1)*drr(1) + drr(2)*drr(2) + drr(3)*drr(3))
ll_eval_energy = ll_eval_energy + K_angle * (1.0 + dr(1)*drr(1) + dr(2)*drr(2) + dr(3)*drr(3))
end do
end do
end function ll_eval_energy
integer function ll_move_atom_1(N, Z, pos, n_extra_data, extra_data, cell, d_i, d_pos, dEmax, dE)
use example_mat_mod
use example_bead_spring_polymer_params_mod
implicit none
integer :: N
integer :: Z(N)
double precision :: pos(3,N), cell(3,3)
integer :: n_extra_data
double precision :: extra_data(n_extra_data, N)
integer :: d_i
double precision :: d_pos(3)
double precision :: dEmax, dE
integer :: i, j
double precision :: dr(3), drp(3), dr_l(3), drp_l(3), dr_l0(3), drp_l0(3), dr_mag, drp_mag, &
dr_mag_sq, drp_mag_sq, pos_l(3,N), d_pos_l(3)
double precision :: cell_inv(3,3)
integer :: dj1, dj2, dj3
double precision, allocatable, save :: new_extra_data(:,:)
integer n_images
double precision cell_height(3), v_norm_hat(3)
print *, "example_bead_spring_polymer_model does not implement move_atom_1 for MC"
call exit(1)
! do i=1, 3
! v_norm_hat = cell(:,mod(i,3)+1) .cross. cell(:,mod(i+1,3)+1)
! v_norm_hat = v_norm_hat / sqrt(sum(v_norm_hat**2))
! cell_height(i) = abs(sum(v_norm_hat*cell(:,i)))
! end do
! n_images = ceiling(cutoff/minval(cell_height))
!
! call matrix3x3_inverse(cell, cell_inv)
! ! into lattice coodinates
! do i=1, N
! pos_l(1,i) = sum(cell_inv(1,:)*pos(:,i))
! pos_l(2,i) = sum(cell_inv(2,:)*pos(:,i))
! pos_l(3,i) = sum(cell_inv(3,:)*pos(:,i))
! end do
! d_pos_l(1) = sum(cell_inv(1,:)*d_pos(:))
! d_pos_l(2) = sum(cell_inv(2,:)*d_pos(:))
! d_pos_l(3) = sum(cell_inv(3,:)*d_pos(:))
!
!
! if (n_extra_data == 1 .and. allocated(new_extra_data)) then
! if (any(shape(new_extra_data) /= shape(extra_data))) then
! deallocate(new_extra_data)
! endif
! endif
! if (n_extra_data == 1 .and. .not. allocated(new_extra_data)) then
! allocate(new_extra_data(n_extra_data, N))
! endif
!
! if (n_extra_data == 1) new_extra_data = extra_data
!
! dE = 0.0
! i=d_i
! do j=1,N
! if (j == i) cycle
!
! dr_l0 = pos_l(:,i) - pos_l(:,j)
! dr_l0 = dr_l0 - floor(dr_l0+0.5)
!
! drp_l0 = pos_l(:,i)+d_pos_l(:) - pos_l(:,j)
! drp_l0 = drp_l0 - floor(drp_l0+0.5)
!
! do dj1=-n_images,n_images
! dr_l(1) = dr_l0(1) + real(dj1, 8)
! drp_l(1) = drp_l0(1) + real(dj1, 8)
! do dj2=-n_images,n_images
! dr_l(2) = dr_l0(2) + real(dj2, 8)
! drp_l(2) = drp_l0(2) + real(dj2, 8)
! do dj3=-n_images,n_images
! dr_l(3) = dr_l0(3) + real(dj3, 8)
! drp_l(3) = drp_l0(3) + real(dj3, 8)
!
! dr(1) = cell(1,1) * dr_l(1) + cell(1,2) * dr_l(2) + cell(1,3) + dr_l(3) ! sum(cell(1,:)*dr_l)
! dr(2) = cell(2,1) * dr_l(1) + cell(2,2) * dr_l(2) + cell(2,3) + dr_l(3) ! sum(cell(2,:)*dr_l)
! dr(3) = cell(3,1) * dr_l(1) + cell(3,2) * dr_l(2) + cell(3,3) + dr_l(3) ! sum(cell(3,:)*dr_l)
!
! drp(1) = cell(1,1) * drp_l(1) + cell(1,2) * drp_l(2) + cell(1,3) + drp_l(3) ! sum(cell(1,:)*drp_l)
! drp(2) = cell(2,1) * drp_l(1) + cell(2,2) * drp_l(2) + cell(2,3) + drp_l(3) ! sum(cell(2,:)*drp_l)
! drp(3) = cell(3,1) * drp_l(1) + cell(3,2) * drp_l(2) + cell(3,3) + drp_l(3) ! sum(cell(3,:)*drp_l)
!
! dr_mag_sq = dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3) ! sum(dr*dr)
! drp_mag_sq = drp(1)*drp(1) + drp(2)*drp(2) + drp(3)*drp(3) ! sum(drp*drp)
!
! if (dr_mag_sq < cutoff_sq) then
! dr_mag = sqrt(dr_mag_sq)
! dE = dE - epsilon*(((sigma/dr_mag)**12 - (sigma/dr_mag)**6) - E_offset)
! if (n_extra_data == 1) then
! new_extra_data(1,i) = new_extra_data(1,i) - 0.5*((sigma/dr_mag)**12 - &
! (sigma/dr_mag)**6 - E_offset)
! new_extra_data(1,j) = new_extra_data(1,j) - 0.5*((sigma/dr_mag)**12 - &
! (sigma/dr_mag)**6 - E_offset)
! endif
! endif
! if (drp_mag_sq < cutoff_sq) then
! drp_mag = sqrt(drp_mag_sq)
! dE = dE + epsilon*(((sigma/drp_mag)**12 - (sigma/drp_mag)**6) - E_offset)
! if (n_extra_data == 1) then
! new_extra_data(1,i) = new_extra_data(1,i) + 0.5*epsilon*(((sigma/drp_mag)**12 - &
! (sigma/drp_mag)**6) - E_offset)
! new_extra_data(1,j) = new_extra_data(1,j) + 0.5*epsilon*(((sigma/drp_mag)**12 - &
! (sigma/drp_mag)**6) - E_offset)
! endif
! endif
!
! end do
! end do
! end do
! end do
!
! if (dE < dEmax) then ! accept
! pos(:,i) = pos(:,i) + d_pos(:)
! if (n_extra_data == 1) extra_data = new_extra_data
! ll_move_atom_1 = 1
! else ! reject
! dE = 0.0
! ll_move_atom_1 = 0
! endif
end function ll_move_atom_1
function ll_eval_forces(N, Z, pos, n_extra_data, extra_data, cell, forces) result(energy)
use example_mat_mod
use example_bead_spring_polymer_params_mod
implicit none
integer :: N
integer :: Z(N)
double precision :: pos(3,N), cell(3,3), forces(3,N)
integer :: n_extra_data
double precision :: extra_data(n_extra_data, N)
double precision :: energy ! result
integer :: i, j
double precision :: dr(3), dr_mag, dr_mag_sq, dr_l(3), dr_l0(3), pos_l(3,N), drr(3), drr_mag
double precision :: cell_inv(3,3), E_term
integer :: dj1, dj2, dj3
integer n_images
double precision cell_height(3), v_norm_hat(3), f_term(3)
integer :: i_chain, j_chain, n_chains, i_bond
do i=1, 3
v_norm_hat = cell(:,mod(i,3)+1) .cross. cell(:,mod(i+1,3)+1)
v_norm_hat = v_norm_hat / sqrt(sum(v_norm_hat**2))
cell_height(i) = abs(sum(v_norm_hat*cell(:,i)))
end do
n_images = ceiling(cutoff/minval(cell_height))
call matrix3x3_inverse(cell, cell_inv)
do i=1, N
pos_l(1,i) = sum(cell_inv(1,:)*pos(:,i))
pos_l(2,i) = sum(cell_inv(2,:)*pos(:,i))
pos_l(3,i) = sum(cell_inv(3,:)*pos(:,i))
end do
if (n_extra_data == 1) extra_data = 0.0
energy = 0.0
forces = 0.0
do i=1, N
i_chain = int((i-1)/chain_length)
do j=i, N
j_chain = int((j-1)/chain_length)
dr_l0 = pos_l(:,i) - pos_l(:,j)
dr_l0 = dr_l0 - floor(dr_l0+0.5)
do dj1=-n_images,n_images
dr_l(1) = dr_l0(1) + real(dj1, 8)
do dj2=-n_images,n_images
dr_l(2) = dr_l0(2) + real(dj2, 8)
do dj3=-n_images,n_images
dr_l(3) = dr_l0(3) + real(dj3, 8)
if (i == j .and. dj1 == 0 .and. dj2 == 0 .and. dj3 == 0) cycle ! no self-interaction
if (i_chain == j_chain .and. abs(i-j) == 1 .and. dj1 == 0 .and. dj2 == 0 .and. dj3 == 0) cycle ! bond excluded
dr(1) = cell(1,1) * dr_l(1) + cell(1,2) * dr_l(2) + cell(1,3) * dr_l(3) ! sum(cell(1,:)*dr_l)
dr(2) = cell(2,1) * dr_l(1) + cell(2,2) * dr_l(2) + cell(2,3) * dr_l(3) ! sum(cell(2,:)*dr_l)
dr(3) = cell(3,1) * dr_l(1) + cell(3,2) * dr_l(2) + cell(3,3) * dr_l(3) ! sum(cell(3,:)*dr_l)
dr_mag_sq = dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3) ! sum(dr*dr)
if (dr_mag_sq < cutoff_sq) then
dr_mag = sqrt(dr_mag_sq)
E_term = epsilon*((sigma/dr_mag)**12 - (sigma/dr_mag)**6 - E_offset)
if (i == j) E_term = E_term * 0.5
energy = energy + E_term
if (n_extra_data == 1) then
extra_data(1,i) = extra_data(1,i) + 0.5*E_term
extra_data(1,j) = extra_data(1,j) + 0.5*E_term
endif
if (i /= j) then
f_term = epsilon*(-12.0*sigma**12/dr_mag**13 + 6.0*sigma**6/dr_mag**7)*(dr/dr_mag)
forces(:,i) = forces(:,i) - f_term
forces(:,j) = forces(:,j) + f_term
endif
endif
end do
end do
end do
end do
end do
n_chains = N/chain_length
do i_chain=1, n_chains
! bonds
do i_bond=(i_chain-1)*chain_length+1, i_chain*chain_length-1
dr = pos(:,i_bond+1) - pos(:, i_bond)
dr_mag_sq = dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3)
dr_mag = sqrt(dr_mag_sq)
energy = energy + K_bond * (dr_mag-r0_bond)*(dr_mag-r0_bond)
f_term = K_bond*2.0*(dr_mag-r0_bond)*dr/dr_mag
forces(:,i_bond) = forces(:,i_bond) + f_term
forces(:,i_bond+1) = forces(:,i_bond+1) - f_term
end do
! angles
do i_bond=(i_chain-1)*chain_length+2, i_chain*chain_length-1
dr = pos(:,i_bond+1) - pos(:, i_bond)
dr_mag = sqrt(dr(1)*dr(1) + dr(2)*dr(2) + dr(3)*dr(3))
dr = dr / dr_mag
drr = pos(:,i_bond-1) - pos(:, i_bond)
drr_mag = sqrt(drr(1)*drr(1) + drr(2)*drr(2) + drr(3)*drr(3))
drr = drr / drr_mag
energy = energy + K_angle * (1.0 + dr(1)*drr(1) + dr(2)*drr(2) + dr(3)*drr(3))
! dr
f_term = K_angle * (drr/dr_mag - (dr(1)*drr(1) + dr(2)*drr(2) + dr(3)*drr(3)) * dr/dr_mag)
forces(:,i_bond) = forces(:, i_bond) + f_term
forces(:,i_bond+1) = forces(:, i_bond+1) - f_term
! drr
f_term = K_angle * (dr/drr_mag - (dr(1)*drr(1) + dr(2)*drr(2) + dr(3)*drr(3)) * drr/drr_mag)
forces(:,i_bond) = forces(:, i_bond) + f_term
forces(:,i_bond-1) = forces(:, i_bond-1) - f_term
end do
end do
end function ll_eval_forces