-
Notifications
You must be signed in to change notification settings - Fork 105
/
Copy pathbd_nvt_lj.f90
276 lines (210 loc) · 12.2 KB
/
bd_nvt_lj.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
! bd_nvt_lj.f90
! Brownian dynamics, NVT ensemble
PROGRAM bd_nvt_lj
!------------------------------------------------------------------------------------------------!
! 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. !
!------------------------------------------------------------------------------------------------!
! Takes in a configuration of atoms (positions, velocities)
! Cubic periodic boundary conditions
! Conducts molecular dynamics using BAOAB algorithm of BJ Leimkuhler and C Matthews
! Appl. Math. Res. eXpress 2013, 34–56 (2013); J. Chem. Phys. 138, 174102 (2013)
! Uses no special neighbour lists
! Reads several variables and options from standard input using a namelist nml
! Leave namelist empty to accept supplied defaults
! Positions r are divided by box length after reading in
! However, input configuration, output configuration, most calculations, and all results
! are given in simulation units defined by the model
! For example, for Lennard-Jones, sigma = 1, epsilon = 1
! We assume mass m=1 throughout
! Despite the program name, there is nothing here specific to Lennard-Jones
! The model is defined in md_module
USE, INTRINSIC :: iso_fortran_env, ONLY : input_unit, output_unit, error_unit, iostat_end, iostat_eor, &
& COMPILER_VERSION, COMPILER_OPTIONS
USE config_io_module, ONLY : read_cnf_atoms, write_cnf_atoms
USE averages_module, ONLY : run_begin, run_end, blk_begin, blk_end, blk_add
USE maths_module, ONLY : random_normals
USE md_module, ONLY : introduction, conclusion, allocate_arrays, deallocate_arrays, &
& force, r, v, f, n, potential_type
IMPLICIT NONE
! Most important variables
REAL :: box ! Box length
REAL :: dt ! Time step
REAL :: r_cut ! Potential cutoff distance
REAL :: temperature ! Temperature (specified)
REAL :: gamma ! Friction coefficient
! Composite interaction = pot & cut & vir & lap & ovr variables
TYPE(potential_type) :: total
INTEGER :: blk, stp, nstep, nblock, ioerr
CHARACTER(len=4), PARAMETER :: cnf_prefix = 'cnf.'
CHARACTER(len=3), PARAMETER :: inp_tag = 'inp'
CHARACTER(len=3), PARAMETER :: out_tag = 'out'
CHARACTER(len=3) :: sav_tag = 'sav' ! May be overwritten with block number
NAMELIST /nml/ nblock, nstep, r_cut, dt, gamma, temperature
WRITE ( unit=output_unit, fmt='(a)' ) 'bd_nvt_lj'
WRITE ( unit=output_unit, fmt='(2a)' ) 'Compiler: ', COMPILER_VERSION()
WRITE ( unit=output_unit, fmt='(2a/)' ) 'Options: ', COMPILER_OPTIONS()
WRITE ( unit=output_unit, fmt='(a)' ) 'Brownian dynamics, constant-NVT ensemble'
WRITE ( unit=output_unit, fmt='(a)' ) 'Particle mass m=1 throughout'
CALL introduction
CALL RANDOM_INIT ( .FALSE., .TRUE. ) ! Initialize random number generator
! Set sensible default run parameters for testing
nblock = 10
nstep = 20000
r_cut = 2.5
dt = 0.005
gamma = 1.0
temperature = 1.0
! Read run parameters from namelist
! Comment out, or replace, this section if you don't like namelists
READ ( unit=input_unit, nml=nml, iostat=ioerr )
IF ( ioerr /= 0 ) THEN
WRITE ( unit=error_unit, fmt='(a,i15)') 'Error reading namelist nml from standard input', ioerr
IF ( ioerr == iostat_eor ) WRITE ( unit=error_unit, fmt='(a)') 'End of record'
IF ( ioerr == iostat_end ) WRITE ( unit=error_unit, fmt='(a)') 'End of file'
STOP 'Error in bd_nvt_lj'
END IF
! Write out run parameters
WRITE ( unit=output_unit, fmt='(a,t40,i15)' ) 'Number of blocks', nblock
WRITE ( unit=output_unit, fmt='(a,t40,i15)' ) 'Number of steps per block', nstep
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Potential cutoff distance', r_cut
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Time step', dt
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Friction coefficient', gamma
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Temperature', temperature
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Ideal diffusion coefft', temperature / gamma
! Read in initial configuration and allocate necessary arrays
CALL read_cnf_atoms ( cnf_prefix//inp_tag, n, box ) ! First call is just to get n and box
WRITE ( unit=output_unit, fmt='(a,t40,i15)' ) 'Number of particles', n
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Simulation box length', box
WRITE ( unit=output_unit, fmt='(a,t40,f15.6)' ) 'Density', REAL(n) / box**3
CALL allocate_arrays ( box, r_cut )
CALL read_cnf_atoms ( cnf_prefix//inp_tag, n, box, r, v ) ! Second call gets r and v
r(:,:) = r(:,:) / box ! Convert positions to box units
r(:,:) = r(:,:) - ANINT ( r(:,:) ) ! Periodic boundaries
! Initial forces, potential, etc plus overlap check
CALL force ( box, r_cut, total )
IF ( total%ovr ) THEN
WRITE ( unit=error_unit, fmt='(a)') 'Overlap in initial configuration'
STOP 'Error in bd_nvt_lj'
END IF
! Initialize arrays for averaging and write column headings
CALL run_begin ( calc_variables() )
DO blk = 1, nblock ! Begin loop over blocks
CALL blk_begin
DO stp = 1, nstep ! Begin loop over steps
CALL b_propagator ( dt/2.0 ) ! B kick half-step
CALL a_propagator ( dt/2.0 ) ! A drift half-step
CALL o_propagator ( dt ) ! O random velocities and friction step
CALL a_propagator ( dt/2.0 ) ! A drift half-step
CALL force ( box, r_cut, total )
IF ( total%ovr ) THEN
WRITE ( unit=error_unit, fmt='(a)') 'Overlap in configuration'
STOP 'Error in bd_nvt_lj'
END IF
CALL b_propagator ( dt/2.0 ) ! B kick half-step
! Calculate and accumulate variables for this step
CALL blk_add ( calc_variables() )
END DO ! End loop over steps
CALL blk_end ( blk ) ! Output block averages
IF ( nblock < 1000 ) WRITE(sav_tag,'(i3.3)') blk ! Number configuration by block
CALL write_cnf_atoms ( cnf_prefix//sav_tag, n, box, r*box, v ) ! Save configuration
END DO ! End loop over blocks
CALL run_end ( calc_variables() ) ! Output run averages
CALL write_cnf_atoms ( cnf_prefix//out_tag, n, box, r*box, v ) ! Write out final configuration
CALL deallocate_arrays
CALL conclusion
CONTAINS
SUBROUTINE a_propagator ( t )
IMPLICIT NONE
REAL, INTENT(in) :: t ! Time over which to propagate (typically dt/2)
! Implements the A propagator (drift)
r(:,:) = r(:,:) + t * v(:,:) / box ! Positions in box=1 units
r(:,:) = r(:,:) - ANINT ( r(:,:) ) ! Periodic boundaries (box=1 units)
END SUBROUTINE a_propagator
SUBROUTINE b_propagator ( t )
IMPLICIT NONE
REAL, INTENT(in) :: t ! Time over which to propagate (typically dt/2)
! Implements the B propagator (kick)
v(:,:) = v(:,:) + t * f(:,:)
END SUBROUTINE b_propagator
SUBROUTINE o_propagator ( t )
USE maths_module, ONLY : expm1
IMPLICIT NONE
REAL, INTENT(in) :: t ! Time over which to propagate (typically dt)
! Implements the O propagator (friction and random contributions)
REAL :: x, c
REAL, DIMENSION(3,n) :: zeta
x = gamma * t
c = -expm1(-2*x) ! 1-exp(-2*x), preserving accuracy for small x
c = SQRT ( c )
CALL random_normals ( 0.0, SQRT(temperature), zeta ) ! Random momenta
v = EXP(-x) * v + c * zeta
END SUBROUTINE o_propagator
FUNCTION calc_variables ( ) RESULT(variables)
USE lrc_module, ONLY : potential_lrc, pressure_lrc
USE averages_module, ONLY : variable_type, msd
IMPLICIT NONE
TYPE(variable_type), DIMENSION(8) :: variables ! The 8 variables listed below
! This function returns all variables of interest in an array, for use in the main program
TYPE(variable_type) :: e_s, p_s, c_s, e_f, p_f, c_f, t_k, t_c
REAL :: kin, fsq, vol, rho
! Preliminary calculations
kin = 0.5*SUM(v**2) ! Total kinetic energy
fsq = SUM(f**2) ! Total squared force
vol = box**3 ! Volume
rho = REAL(n) / vol ! Density
! Variables of interest, of type variable_type, containing three components:
! %val: the instantaneous value
! %nam: used for headings
! %method: indicating averaging method
! If not set below, %method adopts its default value of avg
! The %nam and some other components need only be defined once, at the start of the program,
! but for clarity and readability we assign all the values together below
! Internal energy (cut-and-shifted ) per atom
! Total KE plus cut-and-shifted PE divided by N
e_s = variable_type ( nam = 'E/N cut&shifted', val = (kin+total%pot)/REAL(n) )
! Pressure (cut-and-shifted)
! Ideal gas part plus total virial divided by V
p_s = variable_type ( nam = 'P cut&shifted', val = rho*temperature + total%vir/vol )
! Internal energy (full, including LRC) per atom
! LRC plus total KE plus cut (but not shifted) PE, divided by N
e_f = variable_type ( nam = 'E/N full', val = potential_lrc(rho,r_cut) + (kin+total%cut)/REAL(n) )
! Pressure (full, including LRC)
! LRC + ideal gas part plus total virial divided by V plus LRC
p_f = variable_type ( nam = 'P full', val = pressure_lrc(rho,r_cut) + rho*temperature + total%vir/vol )
! Kinetic temperature
! Momentum is not conserved, hence 3N degrees of freedom
t_k = variable_type ( nam = 'T kinetic', val = 2.0*kin/REAL(3*n) )
! Configurational temperature
! Total squared force divided by total Laplacian
t_c = variable_type ( nam = 'T config', val = fsq/total%lap )
! Heat capacity (cut-and-shifted)
! Total energy divided by temperature and sqrt(N) to make result intensive
c_s = variable_type ( nam = 'Cv/N cut&shifted', val = (kin+total%pot)/(temperature*SQRT(REAL(n))), &
& method = msd, instant = .FALSE. )
! Heat capacity (full)
! Total energy divided by temperature and sqrt(N) to make result intensive; LRC does not contribute
c_f = variable_type ( nam = 'Cv/N full', val = (kin+total%cut)/(temperature*SQRT(REAL(n))), &
& method = msd, instant = .FALSE. )
! Collect together for averaging
variables = [ e_s, p_s, e_f, p_f, t_k, t_c, c_s, c_f ]
END FUNCTION calc_variables
END PROGRAM bd_nvt_lj