-
Notifications
You must be signed in to change notification settings - Fork 18
/
example_Jaglasmooth_model.F90
403 lines (340 loc) · 13 KB
/
example_Jaglasmooth_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
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
! publically accessible things required for interface to pymatnest
!
! subroutine ll_init_model(N_params, params)
! integer :: N_params ! number of parameters
! double precision :: params(N_params) ! list of parameters
!
! initializes potential
!
! subroutine ll_init_config(N, Z, pos, cell, Emax)
! integer :: N ! number of atoms
! integer :: Z(N) ! atomic numbers of atoms
! double precision :: pos(3,N), cell(3,3) ! positions, cell vectors
! double precision :: Emax ! maximum energy for config acceptance
!
! initializes a configuration with energy < Emax
! config will be tested for failure after return
!
! double precision function ll_eval_energy(N, Z, pos, n_extra_data, extra_data, cell)
! integer :: N ! number of atoms
! double precision :: pos(3,N), cell(3,3) ! positions, cell vectors
! integer :: n_extra_data ! width of extra data array
! double precision :: extra_data(n_extra_data, N) ! extra data on output
!
! evaluates energy of a config, sets extra_data, returns energy
!
! integer function ll_move_atom_1(N, pos, n_extra_data, extra_data, cell, d_i, d_pos, dEmax, dE)
! integer :: N ! number of atoms
! double precision :: pos(3,N), cell(3,3) ! positions, cell vectors, on output updated (pos only) to be consistent with acceptance/rejection
! integer :: n_extra_data ! width of extra data array
! double precision :: extra_data(n_extra_data, N) ! extra data on input, on output updated to be consistent with acceptance/rejection
! integer :: d_i ! index of atom to be perturbed, 1-based (called from fortran_MC())
! double precision :: d_pos(3) ! displacement of perturbed atom
! double precision :: dEmax ! maximum change in energy for move acceptance
! double precision :: dE ! on output actual change in energy, 0.0 if move is rejected
!
! moves an atom if dE < dEmax
! if move is accepted, updates pos, extra_data, sets dE
! if move is rejected, nothing is updated, dE set to 0.0
! returns 1 for accept, 0 for reject
!
! double precision function ll_eval_forces(N, pos, n_extra_data, extra_data, cell, forces)
! integer :: N ! number of atoms
! double precision :: pos(3,N), cell(3,3), forces(3,N) ! positions, cell vectors, forces
! integer :: n_extra_data ! width of extra data array
! double precision :: extra_data(n_extra_data, N) ! extra data on output
!
! evaluates forces, sets extra_data
! returns energy
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Jagla potential
subroutine ll_init_model(N_params, params)
use example_Jagla_params_mod
implicit none
integer :: N_params
double precision :: params(N_params)
if (N_params==1) then
write(*,*) "No Jagla parameters found, the default will be used:"
sigma(1) = 1.0d0
d_a(1) = 1.72d0
w_a(1) = -1.0d0
w_r(1) = 3.5d0
cutoff(1)= 3.0d0
cutoff_sq = cutoff*cutoff
else
write(*,*) "The following Jagla parameters will be used:"
sigma(1) = 1.0d0
d_a(1) = params(2) !1.72d0
w_a(1) = -1.0d0
w_r(1) = params(1) !3.5d0
cutoff(1)= params(3) !3.0d0
cutoff_sq = cutoff*cutoff
endif
write(*,*) "Hard sphere",sigma(1),"minimum energy",w_a(1)
write(*,*) "minimum location",d_a(1),"cutoff",cutoff(1)
write(*,*) "repulsive max", w_r(1)
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_Jagla_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)
double precision :: cell_inv(3,3), E_term
integer :: dj1, dj2, dj3
integer :: n_images
double precision cell_height(3), v_norm_hat(3)
call matrix3x3_inverse(cell, cell_inv)
! into lattice coodinates
pos_l = matmul(cell_inv, pos)
if (n_extra_data == 1) extra_data = 0.0d0
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(maxval(cutoff)/minval(cell_height))
ll_eval_energy = 0.0d0
do i=1, N
do j=i, N
dr_l0 = pos_l(:,i)-pos_l(:,j)
dr_l0 = dr_l0 - floor(dr_l0+0.5d0)
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
dr(1) = sum(cell(1,:)*dr_l)
dr(2) = sum(cell(2,:)*dr_l)
dr(3) = sum(cell(3,:)*dr_l)
dr_mag_sq = sum(dr*dr)
if (dr_mag_sq < cutoff_sq(1)) then
dr_mag = sqrt(dr_mag_sq)
if (dr_mag >= d_a(1)) then
E_term = 3*((dr_mag-d_a(1))/(cutoff(1)-d_a(1)))**2-2*((dr_mag-d_a(1))/(cutoff(1)-d_a(1)))**3+w_a(1)
elseif (dr_mag >= sigma(1) ) then
E_term = (3*((d_a(1)-dr_mag)/((d_a(1)-sigma(1))*2))**2-2*((d_a(1)-dr_mag)/((d_a(1)-sigma(1))*2))**3)* &
& ((w_r(1)-w_a(1))*2)+w_a(1)
else
E_term = huge(1.0d0)
endif
if (i == j) E_term = E_term * 0.5d0
ll_eval_energy = ll_eval_energy + E_term
if (n_extra_data == 1) then
extra_data(1,i) = extra_data(1,i) + 0.5d0*E_term
extra_data(1,j) = extra_data(1,j) + 0.5d0*E_term
endif
endif
end do
end do
end do
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_Jagla_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, Z_i, Z_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)
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(maxval(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) then
write(*,*) "ll_move_atom_1 is used, n_extra_data=", n_extra_data
stop
endif
if (n_extra_data == 1) new_extra_data = extra_data
dE = 0.0d0
i=d_i
atoms_c: do j=1,N
if (j == i) cycle
dr_l0 = pos_l(:,i) - pos_l(:,j)
dr_l0 = dr_l0 - floor(dr_l0+0.5d0)
drp_l0 = pos_l(:,i)+d_pos_l(:) - pos_l(:,j)
drp_l0 = drp_l0 - floor(drp_l0+0.5d0)
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) = sum(cell(1,:)*dr_l)
dr(2) = sum(cell(2,:)*dr_l)
dr(3) = sum(cell(3,:)*dr_l)
drp(1) = sum(cell(1,:)*drp_l)
drp(2) = sum(cell(2,:)*drp_l)
drp(3) = sum(cell(3,:)*drp_l)
dr_mag_sq = sum(dr*dr)
drp_mag_sq = sum(drp*drp)
if (dr_mag_sq < cutoff_sq(1)) then
dr_mag = sqrt(dr_mag_sq)
if (dr_mag <= sigma(1)) then
dE = dEmax+1.0
exit atoms_c
endif
if (dr_mag >= d_a(1)) then
dE = dE - (3*((dr_mag-d_a(1))/(cutoff(1)-d_a(1)))**2-2*((dr_mag-d_a(1))/(cutoff(1)-d_a(1)))**3+w_a(1))
else
dE = dE - ((3*((d_a(1)-dr_mag)/((d_a(1)-sigma(1))*2))**2-2*((d_a(1)-dr_mag)/((d_a(1)-sigma(1))*2))**3)* &
& ((w_r(1)-w_a(1))*2)+w_a(1))
endif
endif
if (drp_mag_sq < cutoff_sq(1)) then
drp_mag = sqrt(drp_mag_sq)
if (drp_mag <= sigma(1)) then
dE = dEmax+1.0
exit atoms_c
endif
if (drp_mag >= d_a(1)) then
dE = dE + (3*((drp_mag-d_a(1))/(cutoff(1)-d_a(1)))**2-2*((drp_mag-d_a(1))/(cutoff(1)-d_a(1)))**3+w_a(1))
else
dE = dE + ((3*((d_a(1)-drp_mag)/((d_a(1)-sigma(1))*2))**2-2*((d_a(1)-drp_mag)/((d_a(1)-sigma(1))*2))**3)* &
& ((w_r(1)-w_a(1))*2)+w_a(1))
endif
! if (n_extra_data == 1) then
! new_extra_data(1,i) = new_extra_data(1,i) + 0.5*epsilon(Z_i,Z_j)*(((sigma(Z_i,Z_j)/drp_mag)**12 - &
! (sigma(Z_i,Z_j)/drp_mag)**6) - E_offset(Z_i,Z_j))
! new_extra_data(1,j) = new_extra_data(1,j) + 0.5*epsilon(Z_i,Z_j)*(((sigma(Z_i,Z_j)/drp_mag)**12 - &
! (sigma(Z_i,Z_j)/drp_mag)**6) - E_offset(Z_i,Z_j))
! endif
endif
end do
end do
end do
end do atoms_c
!write(*,*) "move", dE, dEmax
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_LJ_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, Z_i, Z_j
double precision :: dr(3), dr_mag, dr_mag_sq, dr_l(3), dr_l0(3), pos_l(3,N)
double precision :: cell_inv(3,3), E_term
integer :: dj1, dj2, dj3
integer n_images
double precision cell_height(3), v_norm_hat(3)
write(*,*) "Forces are not available for the Jagla potential. Use MC. Stop."
stop
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(maxval(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
Z_i = Z(i)
do j=i, N
Z_j = Z(j)
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
dr(1) = sum(cell(1,:)*dr_l)
dr(2) = sum(cell(2,:)*dr_l)
dr(3) = sum(cell(3,:)*dr_l)
dr_mag_sq = sum(dr*dr)
if (dr_mag_sq < cutoff_sq(Z_i,Z_j)) then
dr_mag = sqrt(dr_mag_sq)
E_term = epsilon(Z_i,Z_j)*((sigma(Z_i,Z_j)/dr_mag)**12 - (sigma(Z_i,Z_j)/dr_mag)**6 - E_offset(Z_i,Z_j))
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
forces(:,i) = forces(:,i) - epsilon(Z_i,Z_j)*(-12.0*sigma(Z_i,Z_j)**12/dr_mag**13 + &
6.0*sigma(Z_i,Z_j)**6/dr_mag**7)*(dr/dr_mag)
forces(:,j) = forces(:,j) + epsilon(Z_i,Z_j)*(-12.0*sigma(Z_i,Z_j)**12/dr_mag**13 + &
6.0*sigma(Z_i,Z_j)**6/dr_mag**7)*(dr/dr_mag)
endif
endif
end do
end do
end do
end do
end do
end function ll_eval_forces