-
Notifications
You must be signed in to change notification settings - Fork 103
/
mc_lj_module.f90
390 lines (294 loc) · 14.2 KB
/
mc_lj_module.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
! mc_lj_module.f90
! Energy and move routines for MC simulation, LJ potential
MODULE mc_module
!------------------------------------------------------------------------------------------------!
! This software was written in 2016/17 !
! by Michael P. Allen <[email protected]>/<[email protected]> !
! and Dominic J. Tildesley <[email protected]> ("the authors"), !
! to accompany the book "Computer Simulation of Liquids", second edition, 2017 ("the text"), !
! published by Oxford University Press ("the publishers"). !
! !
! LICENCE !
! Creative Commons CC0 Public Domain Dedication. !
! To the extent possible under law, the authors have dedicated all copyright and related !
! and neighboring rights to this software to the PUBLIC domain worldwide. !
! This software is distributed without any warranty. !
! You should have received a copy of the CC0 Public Domain Dedication along with this software. !
! If not, see <http://creativecommons.org/publicdomain/zero/1.0/>. !
! !
! DISCLAIMER !
! The authors and publishers make no warranties about the software, and disclaim liability !
! for all uses of the software, to the fullest extent permitted by applicable law. !
! The authors and publishers do not recommend use of this software for any purpose. !
! It is made freely available, solely to clarify points made in the text. When using or citing !
! the software, you should not imply endorsement by the authors or publishers. !
!------------------------------------------------------------------------------------------------!
USE, INTRINSIC :: iso_fortran_env, ONLY : output_unit, error_unit
IMPLICIT NONE
PRIVATE
! Public routines
PUBLIC :: introduction, conclusion, allocate_arrays, extend_arrays, deallocate_arrays
PUBLIC :: box_check, potential_1, potential, force_sq
PUBLIC :: move, create, destroy
! Public data
INTEGER, PUBLIC :: n ! Number of atoms
REAL, DIMENSION(:,:), ALLOCATABLE, PUBLIC :: r ! Positions (3,n)
! Private data
REAL, DIMENSION(:,:), ALLOCATABLE :: f ! Forces for force_sq calculation (3,n)
REAL, DIMENSION(:,:), ALLOCATABLE :: tmp_r ! Temporary array for dynamic reallocation of r
INTEGER, PARAMETER :: lt = -1, gt = 1 ! Options for j-range
! Public derived type
TYPE, PUBLIC :: potential_type ! A composite variable for interactions comprising
REAL :: pot ! the potential energy cut at r_cut and
REAL :: vir ! the virial and
REAL :: lap ! the Laplacian and
LOGICAL :: ovr ! a flag indicating overlap (i.e. pot too high to use)
CONTAINS
PROCEDURE :: add_potential_type
PROCEDURE :: subtract_potential_type
GENERIC :: OPERATOR(+) => add_potential_type
GENERIC :: OPERATOR(-) => subtract_potential_type
END TYPE potential_type
CONTAINS
FUNCTION add_potential_type ( a, b ) RESULT (c)
IMPLICIT NONE
TYPE(potential_type) :: c ! Result is the sum of
CLASS(potential_type), INTENT(in) :: a, b ! the two inputs
c%pot = a%pot + b%pot
c%vir = a%vir + b%vir
c%lap = a%lap + b%lap
c%ovr = a%ovr .OR. b%ovr
END FUNCTION add_potential_type
FUNCTION subtract_potential_type ( a, b ) RESULT (c)
IMPLICIT NONE
TYPE(potential_type) :: c ! Result is the difference of
CLASS(potential_type), INTENT(in) :: a, b ! the two inputs
c%pot = a%pot - b%pot
c%vir = a%vir - b%vir
c%lap = a%lap - b%lap
c%ovr = a%ovr .OR. b%ovr ! This is meaningless, but inconsequential
END FUNCTION subtract_potential_type
SUBROUTINE introduction
IMPLICIT NONE
WRITE ( unit=output_unit, fmt='(a)' ) 'Lennard-Jones potential'
WRITE ( unit=output_unit, fmt='(a)' ) 'Cut (but not shifted)'
WRITE ( unit=output_unit, fmt='(a)' ) 'Diameter, sigma = 1'
WRITE ( unit=output_unit, fmt='(a)' ) 'Well depth, epsilon = 1'
END SUBROUTINE introduction
SUBROUTINE conclusion
IMPLICIT NONE
WRITE ( unit=output_unit, fmt='(a)') 'Program ends'
END SUBROUTINE conclusion
SUBROUTINE allocate_arrays ( box, r_cut )
IMPLICIT NONE
REAL, INTENT(in) :: box ! Simulation box length
REAL, INTENT(in) :: r_cut ! Potential cutoff distance
REAL :: r_cut_box
ALLOCATE ( r(3,n), f(3,n) )
r_cut_box = r_cut / box
IF ( r_cut_box > 0.5 ) THEN
WRITE ( unit=error_unit, fmt='(a,f15.6)') 'r_cut/box too large ', r_cut_box
STOP 'Error in allocate_arrays'
END IF
END SUBROUTINE allocate_arrays
SUBROUTINE extend_arrays ( n_new )
IMPLICIT NONE
INTEGER, INTENT(in) :: n_new ! New size of arrays
INTEGER :: n_old
n_old = SIZE ( r, dim=2 )
IF ( n_new <= n_old ) THEN
WRITE ( unit=error_unit, fmt='(a,2i15)') 'Unexpected array size ', n_old, n_new
STOP 'Error in extend_arrays'
END IF
! Reallocate r array, preserving existing data
ALLOCATE ( tmp_r(3,n_new) )
tmp_r(:,1:n_old) = r
tmp_r(:,n_old+1:) = 0.0
CALL MOVE_ALLOC ( tmp_r, r )
! Reallocate f array, no need to preserve data
DEALLOCATE ( f )
ALLOCATE ( f(3,n_new) )
END SUBROUTINE extend_arrays
SUBROUTINE deallocate_arrays
IMPLICIT NONE
DEALLOCATE ( r, f )
END SUBROUTINE deallocate_arrays
SUBROUTINE box_check ( box, r_cut )
IMPLICIT NONE
REAL, INTENT(in) :: box ! Simulation box length
REAL, INTENT(in) :: r_cut ! Potential cutoff distance
REAL :: r_cut_box
r_cut_box = r_cut / box
IF ( r_cut_box > 0.5 ) THEN
WRITE ( unit=error_unit, fmt='(a,f15.6)') 'r_cut/box too large ', r_cut_box
STOP 'Error in box_check'
END IF
END SUBROUTINE box_check
FUNCTION potential ( box, r_cut ) RESULT ( total )
IMPLICIT NONE
TYPE(potential_type) :: total ! Returns a composite of pot, vir etc
REAL, INTENT(in) :: box ! Simulation box length
REAL, INTENT(in) :: r_cut ! Potential cutoff distance
! total%pot is the nonbonded cut (not shifted) potential energy for whole system
! total%vir is the corresponding virial for whole system
! total%lap is the corresponding Laplacian for whole system
! total%ovr is a flag indicating overlap (potential too high) to avoid overflow
! If this flag is .true., the values of total%pot etc should not be used
! Actual calculation is performed by function potential_1
! We assume that the box length may vary, so we check r_cut/box at the start.
TYPE(potential_type) :: partial ! Atomic contribution to total
INTEGER :: i
REAL :: r_cut_box
IF ( n > SIZE(r,dim=2) ) THEN ! should never happen
WRITE ( unit=error_unit, fmt='(a,2i15)' ) 'Array bounds error for r', n, SIZE(r,dim=2)
STOP 'Error in potential'
END IF
r_cut_box = r_cut / box
IF ( r_cut_box > 0.5 ) THEN
WRITE ( unit=error_unit, fmt='(a,f15.6)') 'r_cut/box too large ', r_cut_box
STOP 'Error in potential'
END IF
total = potential_type ( pot=0.0, vir=0.0, lap=0.0, ovr=.FALSE. ) ! Initialize
DO i = 1, n - 1
partial = potential_1 ( r(:,i), i, box, r_cut, gt )
IF ( partial%ovr ) THEN
total%ovr = .TRUE. ! Overlap detected
RETURN ! Return immediately
END IF
total = total + partial
END DO
total%ovr = .FALSE. ! No overlaps detected (redundant, but for clarity)
END FUNCTION potential
FUNCTION potential_1 ( ri, i, box, r_cut, j_range ) RESULT ( partial )
IMPLICIT NONE
TYPE(potential_type) :: partial ! Returns a composite of pot, vir etc for given atom
REAL, DIMENSION(3), INTENT(in) :: ri ! Coordinates of atom of interest
INTEGER, INTENT(in) :: i ! Index of atom of interest
REAL, INTENT(in) :: box ! Simulation box length
REAL, INTENT(in) :: r_cut ! Potential cutoff distance
INTEGER, OPTIONAL, INTENT(in) :: j_range ! Optional partner index range
! partial%pot is the nonbonded cut (not shifted) potential energy of atom ri with a set of other atoms
! partial%vir is the corresponding virial of atom ri
! partial%lap is the corresponding Laplacian of atom ri
! partial%ovr is a flag indicating overlap (potential too high) to avoid overflow
! If this is .true., the values of partial%pot etc should not be used
! The coordinates in ri are not necessarily identical with those in r(:,i)
! The optional argument j_range restricts partner indices to j>i, or j<i
! It is assumed that r has been divided by box
! Results are in LJ units where sigma = 1, epsilon = 1
! We assume that the box length may vary, but that r_cut/box
! has been checked since the last change in box length.
INTEGER :: j, j1, j2
REAL :: r_cut_box, r_cut_box_sq, box_sq
REAL :: sr2, sr6, sr12, rij_sq
REAL, DIMENSION(3) :: rij
REAL, PARAMETER :: sr2_ovr = 1.77 ! overlap threshold (pot > 100)
TYPE(potential_type) :: pair
IF ( n > SIZE(r,dim=2) ) THEN ! should never happen
WRITE ( unit=error_unit, fmt='(a,2i15)' ) 'Array bounds error for r', n, SIZE(r,dim=2)
STOP 'Error in potential_1'
END IF
IF ( PRESENT ( j_range ) ) THEN
SELECT CASE ( j_range )
CASE ( lt ) ! j < i
j1 = 1
j2 = i-1
CASE ( gt ) ! j > i
j1 = i+1
j2 = n
CASE default ! should never happen
WRITE ( unit = error_unit, fmt='(a,i10)') 'j_range error ', j_range
STOP 'Impossible error in potential_1'
END SELECT
ELSE
j1 = 1
j2 = n
END IF
r_cut_box = r_cut / box
r_cut_box_sq = r_cut_box**2
box_sq = box**2
partial = potential_type ( pot=0.0, vir=0.0, lap=0.0, ovr=.FALSE. ) ! Initialize
DO j = j1, j2 ! Loop over selected range of partners
IF ( i == j ) CYCLE ! Skip self
rij(:) = ri(:) - r(:,j) ! Separation vector
rij(:) = rij(:) - ANINT ( rij(:) ) ! Periodic boundaries in box=1 units
rij_sq = SUM ( rij**2 ) ! Squared separation in box=1 units
IF ( rij_sq < r_cut_box_sq ) THEN ! Check within range
rij_sq = rij_sq * box_sq ! Now in sigma=1 units
sr2 = 1.0 / rij_sq ! (sigma/rij)**2
pair%ovr = sr2 > sr2_ovr ! Overlap if too close
IF ( pair%ovr ) THEN
partial%ovr = .TRUE. ! Overlap detected
RETURN ! Return immediately
END IF
sr6 = sr2**3
sr12 = sr6**2
pair%pot = sr12 - sr6 ! LJ pair potential (cut but not shifted)
pair%vir = pair%pot + sr12 ! LJ pair virial
pair%lap = ( 22.0*sr12 - 5.0*sr6 ) * sr2 ! LJ pair Laplacian
partial = partial + pair
END IF ! End check within range
END DO ! End loop over selected range of partners
! Include numerical factors
partial%pot = partial%pot * 4.0 ! 4*epsilon
partial%vir = partial%vir * 24.0 / 3.0 ! 24*epsilon and divide virial by 3
partial%lap = partial%lap * 24.0 * 2.0 ! 24*epsilon and factor 2 for ij and ji
partial%ovr = .FALSE. ! No overlaps detected (redundant, but for clarity)
END FUNCTION potential_1
FUNCTION force_sq ( box, r_cut ) RESULT ( fsq )
IMPLICIT NONE
REAL :: fsq ! Returns total squared force
REAL, INTENT(in) :: box ! Simulation box length
REAL, INTENT(in) :: r_cut ! Potential cutoff distance
! Calculates total squared force (using array f)
INTEGER :: i, j
REAL :: r_cut_box, r_cut_box_sq, box_sq, rij_sq
REAL :: sr2, sr6, sr12
REAL, DIMENSION(3) :: rij, fij
r_cut_box = r_cut / box
r_cut_box_sq = r_cut_box ** 2
box_sq = box ** 2
IF ( r_cut_box > 0.5 ) THEN
WRITE ( unit=error_unit, fmt='(a,f15.6)') 'r_cut/box too large ', r_cut_box
STOP 'Error in force_sq'
END IF
f = 0.0 ! Initialize
DO i = 1, n - 1 ! Begin outer loop over atoms
DO j = i + 1, n ! Begin inner loop over atoms
rij(:) = r(:,i) - r(:,j) ! Separation vector
rij(:) = rij(:) - ANINT ( rij(:) ) ! Periodic boundary conditions in box=1 units
rij_sq = SUM ( rij**2 ) ! Squared separation
IF ( rij_sq < r_cut_box_sq ) THEN ! Check within cutoff
rij_sq = rij_sq * box_sq ! Now in sigma=1 units
rij(:) = rij(:) * box ! Now in sigma=1 units
sr2 = 1.0 / rij_sq
sr6 = sr2 ** 3
sr12 = sr6 ** 2
fij = rij * (2.0*sr12 - sr6) *sr2 ! LJ pair forces
f(:,i) = f(:,i) + fij
f(:,j) = f(:,j) - fij
END IF ! End check within cutoff
END DO ! End inner loop over atoms
END DO ! End outer loop over atoms
f = f * 24.0 ! Numerical factor 24*epsilon
fsq = SUM ( f**2 ) ! Result
END FUNCTION force_sq
SUBROUTINE move ( i, ri )
IMPLICIT NONE
INTEGER, INTENT(in) :: i
REAL, DIMENSION(3), INTENT(in) :: ri
r(:,i) = ri ! New position
END SUBROUTINE move
SUBROUTINE create ( ri )
IMPLICIT NONE
REAL, DIMENSION(3), INTENT(in) :: ri
n = n+1 ! Increase number of atoms
r(:,n) = ri ! Add new atom at the end
END SUBROUTINE create
SUBROUTINE destroy ( i )
IMPLICIT NONE
INTEGER, INTENT(in) :: i
r(:,i) = r(:,n) ! Replace atom i coordinates with atom n
n = n - 1 ! Reduce number of atoms
END SUBROUTINE destroy
END MODULE mc_module