| 1 |
module long_range_forces |
| 2 |
use definitions, ONLY : dp, ndim |
| 3 |
use simulation |
| 4 |
implicit none |
| 5 |
|
| 6 |
|
| 7 |
|
| 8 |
|
| 9 |
real(kind = dp ) :: rlistsq |
| 10 |
real(kind = dp ) :: rcutsq |
| 11 |
|
| 12 |
|
| 13 |
|
| 14 |
|
| 15 |
integer, allocatable, dimension(:) :: point |
| 16 |
integer, allocatable, dimension(:) :: list |
| 17 |
|
| 18 |
|
| 19 |
|
| 20 |
|
| 21 |
contains |
| 22 |
|
| 23 |
|
| 24 |
|
| 25 |
|
| 26 |
subroutine long_range_force ( check_exclude, & |
| 27 |
pair_i, pair_j, n_pairs, natoms, sigma, epslon, & |
| 28 |
box, & |
| 29 |
rx, ry, rz, fx, fy, fz, & |
| 30 |
v, & |
| 31 |
update, point, list, maxnab, & |
| 32 |
rx0, ry0, rz0, & |
| 33 |
ux, uy, uz, & |
| 34 |
tx, ty, tz, & |
| 35 |
ex, ey, ez, & |
| 36 |
mu, rrf, rtaper, & |
| 37 |
dielectric, & |
| 38 |
is_dipole, is_vdw, is_lj, is_ssd, & |
| 39 |
Axx, Axy, Axz, & |
| 40 |
Ayx, Ayy, Ayz, & |
| 41 |
Azx, Azy, Azz ) |
| 42 |
|
| 43 |
implicit none |
| 44 |
|
| 45 |
! This routine implements the Verlet neighbor list for all of |
| 46 |
! the long range force calculations. After a pair has been |
| 47 |
! determined to lie within the neighbor list, the pair is passed |
| 48 |
! to the appropriate force functions. |
| 49 |
|
| 50 |
|
| 51 |
|
| 52 |
! Passed Arguments |
| 53 |
|
| 54 |
integer :: n_pairs ! the number of excluded pairs |
| 55 |
integer :: natoms ! the number of atoms in the system |
| 56 |
integer :: maxnab ! the max number of neighbors for the verlet list |
| 57 |
|
| 58 |
logical(kind=2) :: check_exclude ! boolean to check for exclusion pairs |
| 59 |
logical(kind=2) :: update ! boolean to update the neighbor list |
| 60 |
|
| 61 |
double precision rcut ! the VDW/LJ cutoff radius |
| 62 |
double precision rlist ! Verlet list cutoff |
| 63 |
double precision box_x, box_y, box_z ! periodic boundry conditions |
| 64 |
double precision rrf ! dipole reaction field cutoff |
| 65 |
double precision rtaper ! the taper radius of the rxn field |
| 66 |
double precision dielectric ! the dielectric cnst. of the rxn field |
| 67 |
double precision v ! the potential energy |
| 68 |
|
| 69 |
|
| 70 |
|
| 71 |
! Passed Arrays |
| 72 |
|
| 73 |
integer, dimension(n_pairs):: pair_i, pair_j ! the excluded pairs i and j |
| 74 |
integer, dimension(natoms) :: point ! the Verlet list ptr array |
| 75 |
integer, dimension(maxnab) :: list ! the verlet list |
| 76 |
|
| 77 |
logical(kind=2), dimension(natoms) :: is_dipole ! dipole boolean array |
| 78 |
logical(kind=2), dimension(natoms) :: is_vdw ! VDW boolean array |
| 79 |
logical(kind=2), dimension(natoms) :: is_lj ! LJ boolean array |
| 80 |
logical(kind=2), dimension(natoms) :: is_ssd ! ssd boolean array |
| 81 |
|
| 82 |
double precision, dimension(natoms) :: sigma ! VDW/LJ distance prmtr. |
| 83 |
double precision, dimension(natoms) :: epslon ! VDW/LJ well depth prmtr. |
| 84 |
double precision, dimension(natoms) :: rx, ry, rz ! positions |
| 85 |
double precision, dimension(natoms) :: fx, fy, fz ! forces |
| 86 |
double precision, dimension(natoms) :: rx0, ry0, rz0 ! intial verlet positions |
| 87 |
double precision, dimension(natoms) :: ux, uy, uz ! dipole unit vectors |
| 88 |
double precision, dimension(natoms) :: tx, ty, tz ! torsions |
| 89 |
double precision, dimension(natoms) :: ex, ey, ez ! reacion field |
| 90 |
double precision, dimension(natoms) :: mu ! dipole moments |
| 91 |
double precision, dimension(natoms) :: Axx, Axy, Axz ! rotational array |
| 92 |
double precision, dimension(natoms) :: Ayx, Ayy, Ayz ! |
| 93 |
double precision, dimension(natoms) :: Azx, Azy, Azz ! |
| 94 |
|
| 95 |
! local variables |
| 96 |
|
| 97 |
double precision rxi, ryi, rzi |
| 98 |
double precision rcutsq, rlstsq, rrfsq, rijsq, rij |
| 99 |
double precision rxij, ryij, rzij |
| 100 |
double precision prerf ! a reaction field preterm |
| 101 |
|
| 102 |
double precision, dimension(9,natoms) :: A |
| 103 |
|
| 104 |
integer :: i, j, k, ix, jx, nlist |
| 105 |
integer :: jbeg, jend, jnab |
| 106 |
|
| 107 |
logical :: exclude_temp1, exclude_temp2, exclude |
| 108 |
|
| 109 |
|
| 110 |
!******************************************************************* |
| 111 |
|
| 112 |
do i=1, natoms |
| 113 |
A(1,i) = Axx(i) |
| 114 |
A(2,i) = Axy(i) |
| 115 |
A(3,i) = Axz(i) |
| 116 |
|
| 117 |
A(4,i) = Ayx(i) |
| 118 |
A(5,i) = Ayy(i) |
| 119 |
A(6,i) = Ayz(i) |
| 120 |
|
| 121 |
A(7,i) = Azx(i) |
| 122 |
A(8,i) = Azy(i) |
| 123 |
A(9,i) = Azz(i) |
| 124 |
end do |
| 125 |
|
| 126 |
|
| 127 |
rlstsq = rlist * rlist |
| 128 |
rcutsq = rcut * rcut |
| 129 |
rrfsq = rrf * rrf |
| 130 |
|
| 131 |
prerf = 14.39257d0 * 2.0d0 * ( dielectric - 1.0d0 ) / & |
| 132 |
( ( 2.0d0 * dielectric + 1.0d0 ) * rrfsq * rrf ) |
| 133 |
|
| 134 |
! ** zero forces ** |
| 135 |
|
| 136 |
do I = 1, natoms |
| 137 |
|
| 138 |
fx(i) = 0.0d0 |
| 139 |
fy(i) = 0.0d0 |
| 140 |
fz(i) = 0.0d0 |
| 141 |
|
| 142 |
tx(i) = 0.0d0 |
| 143 |
ty(i) = 0.0d0 |
| 144 |
tz(i) = 0.0d0 |
| 145 |
|
| 146 |
ex(i) = 0.0d0 |
| 147 |
ey(i) = 0.0d0 |
| 148 |
ez(i) = 0.0d0 |
| 149 |
|
| 150 |
end do |
| 151 |
|
| 152 |
V = 0.0d0 |
| 153 |
|
| 154 |
if ( update ) then |
| 155 |
|
| 156 |
|
| 157 |
! ** save current configuration, construct ** |
| 158 |
! ** neighbour list and calculate forces ** |
| 159 |
|
| 160 |
call long_range_save(natoms,rx,ry,rz,rx0,ry0,rz0) |
| 161 |
|
| 162 |
nlist = 0 |
| 163 |
|
| 164 |
do i = 1, natoms - 1 |
| 165 |
|
| 166 |
point(i) = nlist + 1 |
| 167 |
|
| 168 |
rxi = rx(i) |
| 169 |
ryi = ry(i) |
| 170 |
rzi = rz(i) |
| 171 |
|
| 172 |
do j = i + 1, natoms |
| 173 |
|
| 174 |
rxij = rx(j) - rxi |
| 175 |
ryij = ry(j) - ryi |
| 176 |
rzij = rz(j) - rzi |
| 177 |
|
| 178 |
exclude = .false. |
| 179 |
|
| 180 |
if( check_exclude ) then |
| 181 |
|
| 182 |
do k = 1, n_pairs |
| 183 |
|
| 184 |
! add 1 for c -> fortran indexing |
| 185 |
ix = pair_i(k) + 1 |
| 186 |
jx = pair_j(k) + 1 |
| 187 |
exclude_temp1 = (i.eq.ix).and.(j.eq.jx) |
| 188 |
exclude_temp2 = (i.eq.jx).and.(j.eq.ix) |
| 189 |
|
| 190 |
if(exclude_temp1.or.exclude_temp2) then |
| 191 |
|
| 192 |
exclude = .true. |
| 193 |
endif |
| 194 |
enddo |
| 195 |
endif |
| 196 |
|
| 197 |
if(.not.exclude) then |
| 198 |
|
| 199 |
rxij = rxij - box_x * dsign( 1.0d0, rxij ) & |
| 200 |
* int( (dabs( rxij / box_x ) + 0.5d0) ) |
| 201 |
ryij = ryij - box_y * dsign( 1.0d0, ryij ) & |
| 202 |
* int( (dabs( ryij / box_y ) + 0.5d0) ) |
| 203 |
rzij = rzij - box_z * dsign( 1.0d0, rzij ) & |
| 204 |
* int( (dabs( rzij / box_z ) + 0.5d0) ) |
| 205 |
|
| 206 |
rijsq = rxij * rxij + ryij * ryij + rzij * rzij |
| 207 |
|
| 208 |
if ( rijsq .lt. rlstsq ) then |
| 209 |
|
| 210 |
nlist = nlist + 1 |
| 211 |
list(nlist) = j |
| 212 |
|
| 213 |
! ** remove this check if maxnab is appropriate ** |
| 214 |
|
| 215 |
if ( nlist .eq. maxnab ) then |
| 216 |
write(*,*) i, j, nlist, maxnab, rijsq, rlstsq |
| 217 |
stop 'list too small' |
| 218 |
endif |
| 219 |
|
| 220 |
if ( rijsq .lt. rcutsq ) then |
| 221 |
|
| 222 |
if ( is_vdw(i) .and. is_vdw(j) ) then |
| 223 |
call force_VDW( i, j, rcutsq, rijsq, sigma, epslon, & |
| 224 |
v, fx, fy, fz, rxij, ryij, rzij, natoms ) |
| 225 |
endif |
| 226 |
|
| 227 |
if ( is_lj(i) .and. is_lj(j) ) then |
| 228 |
call force_lj( i, j, rcutsq, rijsq, sigma, epslon, v, & |
| 229 |
fx, fy, fz, rxij, ryij, rzij, natoms ) |
| 230 |
endif |
| 231 |
|
| 232 |
if( is_ssd(i) .and. is_ssd(j) ) then |
| 233 |
|
| 234 |
rij = dsqrt(rijsq); |
| 235 |
|
| 236 |
call ssd_forces( natoms, i, j, i, j, rxij, ryij, rzij, & |
| 237 |
rij, v, rijsq, fx, fy, fz, tx, ty, tz, A ) |
| 238 |
end if |
| 239 |
|
| 240 |
endif |
| 241 |
|
| 242 |
if( is_dipole(i) .and. is_dipole(j) ) then |
| 243 |
if( rijsq .lt. rrfsq ) then |
| 244 |
|
| 245 |
call dipole_rf( i, j, rrfsq, rijsq, v, rxij, ryij, & |
| 246 |
rzij, rrf, dielectric, rtaper, natoms, & |
| 247 |
ux, uy, uz, mu, fx, fy, fz, tx, ty, tz, & |
| 248 |
ex, ey, ez, prerf ) |
| 249 |
endif |
| 250 |
end if |
| 251 |
endif |
| 252 |
endif |
| 253 |
|
| 254 |
enddo |
| 255 |
|
| 256 |
end do |
| 257 |
|
| 258 |
point(natoms) = nlist + 1 |
| 259 |
|
| 260 |
else |
| 261 |
|
| 262 |
|
| 263 |
!** use the list to find the neighbours ** |
| 264 |
|
| 265 |
do i = 1, natoms - 1 |
| 266 |
|
| 267 |
jbeg = point(i) |
| 268 |
jend = point(i+1) - 1 |
| 269 |
|
| 270 |
!** check that atom i has neighbours ** |
| 271 |
|
| 272 |
if( jbeg .le. jend ) then |
| 273 |
|
| 274 |
rxi = rx(i) |
| 275 |
ryi = ry(i) |
| 276 |
rzi = rz(i) |
| 277 |
|
| 278 |
do jnab = jbeg, jend |
| 279 |
|
| 280 |
j = list(jnab) |
| 281 |
|
| 282 |
rxij = rx(j) - rxi |
| 283 |
ryij = ry(j) - ryi |
| 284 |
rzij = rz(j) - rzi |
| 285 |
|
| 286 |
rxij = rxij - box_x * dsign( 1.0d0, rxij ) & |
| 287 |
* int( dabs( rxij / box_x ) + 0.5d0 ) |
| 288 |
ryij = ryij - box_y * dsign( 1.0d0, ryij ) & |
| 289 |
* int( dabs( ryij / box_y ) + 0.5d0 ) |
| 290 |
rzij = rzij - box_z * dsign( 1.0d0, rzij ) & |
| 291 |
* int( dabs( rzij / box_z ) + 0.5d0 ) |
| 292 |
|
| 293 |
rijsq = rxij * rxij + ryij * ryij + rzij * rzij |
| 294 |
|
| 295 |
if ( rijsq .lt. rcutsq ) then |
| 296 |
|
| 297 |
if ( is_vdw(i) .and. is_vdw(j) ) then |
| 298 |
call force_VDW( i, j, rcutsq, rijsq, sigma, epslon, & |
| 299 |
v, fx, fy, fz, rxij, ryij, rzij, natoms ) |
| 300 |
endif |
| 301 |
|
| 302 |
if ( is_lj(i) .and. is_lj(j) ) then |
| 303 |
call force_lj( i, j, rcutsq, rijsq, sigma, epslon, v, & |
| 304 |
fx, fy, fz, rxij, ryij, rzij, natoms ) |
| 305 |
endif |
| 306 |
|
| 307 |
if( is_ssd(i) .and. is_ssd(j) ) then |
| 308 |
|
| 309 |
rij = dsqrt(rijsq); |
| 310 |
|
| 311 |
call ssd_forces( natoms, i, j, i, j, rxij, ryij, rzij, & |
| 312 |
rij, v, rijsq, fx, fy, fz, tx, ty, tz, A ) |
| 313 |
end if |
| 314 |
|
| 315 |
endif |
| 316 |
|
| 317 |
|
| 318 |
if( is_dipole(i) .and. is_dipole(j) ) then |
| 319 |
if( rijsq .lt. rrfsq ) then |
| 320 |
|
| 321 |
call dipole_rf( i, j, rrfsq, rijsq, v, rxij, ryij, & |
| 322 |
rzij, rrf, dielectric, rtaper, natoms, & |
| 323 |
ux, uy, uz, mu, fx, fy, fz, tx, ty, tz, & |
| 324 |
ex, ey, ez, prerf ) |
| 325 |
endif |
| 326 |
endif |
| 327 |
|
| 328 |
end do |
| 329 |
endif |
| 330 |
end do |
| 331 |
endif |
| 332 |
|
| 333 |
! calculate the reaction field effect |
| 334 |
|
| 335 |
!!$ call reaction_field( v, prerf, natoms, mu, ux, uy, uz, tx, ty, tz, & |
| 336 |
!!$ ex, ey, ez, is_dipole ) |
| 337 |
|
| 338 |
|
| 339 |
end subroutine long_range_force |
| 340 |
|
| 341 |
|
| 342 |
|
| 343 |
subroutine long_range_check ( rcut, rlist, update, natoms, & |
| 344 |
rx,ry,rz,rx0,ry0,rz0) |
| 345 |
|
| 346 |
|
| 347 |
!******************************************************************* |
| 348 |
!** decides whether the list needs to be reconstructed. ** |
| 349 |
!** ** |
| 350 |
!** principal variables: ** |
| 351 |
!** ** |
| 352 |
!** real rx(n),ry(n),rz(n) atom positions ** |
| 353 |
!** real rx0(n),ry0(n),rz0(n) coordinates at last update ** |
| 354 |
!** real rlist radius of verlet list ** |
| 355 |
!** real rcut cutoff distance for forces ** |
| 356 |
!** real dispmx largest displacement ** |
| 357 |
!** integer n number of atoms ** |
| 358 |
!** logical update if true the list is updated ** |
| 359 |
!** ** |
| 360 |
!** usage: ** |
| 361 |
!** ** |
| 362 |
!** check is called to set update before every call to force. ** |
| 363 |
!******************************************************************* |
| 364 |
|
| 365 |
|
| 366 |
! Passed variables |
| 367 |
|
| 368 |
integer :: natoms |
| 369 |
|
| 370 |
logical(kind=2) :: update |
| 371 |
|
| 372 |
double precision rcut, rlist |
| 373 |
|
| 374 |
! Passed arrays |
| 375 |
double precision, dimension(natoms) :: rx, ry, rz |
| 376 |
double precision, dimension(natoms) :: rx0, ry0, rz0 |
| 377 |
|
| 378 |
! local variables |
| 379 |
|
| 380 |
double precision dispmx |
| 381 |
|
| 382 |
integer :: i |
| 383 |
|
| 384 |
!******************************************************************* |
| 385 |
|
| 386 |
!** calculate maximum displacement since last update ** |
| 387 |
|
| 388 |
dispmx = 0.0d0 |
| 389 |
|
| 390 |
do i = 1, natoms |
| 391 |
|
| 392 |
!write(*,*) 'calling displacement', i |
| 393 |
|
| 394 |
dispmx = max ( dabs ( rx(i) - rx0(i) ), dispmx ) |
| 395 |
dispmx = max ( dabs ( ry(i) - ry0(i) ), dispmx ) |
| 396 |
dispmx = max ( dabs ( rz(i) - rz0(i) ), dispmx ) |
| 397 |
|
| 398 |
!write(*,*) 'called displacement', i |
| 399 |
|
| 400 |
end do |
| 401 |
|
| 402 |
!** a conservative test of the list skin crossing ** |
| 403 |
|
| 404 |
dispmx = 2.0d0 * dsqrt ( 3.0d0 * dispmx ** 2 ) |
| 405 |
|
| 406 |
update = ( dispmx .gt. ( rlist - rcut ) ) |
| 407 |
|
| 408 |
|
| 409 |
return |
| 410 |
end subroutine long_range_check |
| 411 |
|
| 412 |
|
| 413 |
|
| 414 |
subroutine long_range_save(natoms,rx,ry,rz,rx0,ry0,rz0) |
| 415 |
|
| 416 |
|
| 417 |
!******************************************************************* |
| 418 |
!** saves current configuration for future checking. ** |
| 419 |
!** ** |
| 420 |
!** principal variables: ** |
| 421 |
!** ** |
| 422 |
!** real rx(n),ry(n),rz(n) atom positions ** |
| 423 |
!** real rx0(n),ry0(n),rz0(n) storage locations for save ** |
| 424 |
!** integer n number of atoms ** |
| 425 |
!** ** |
| 426 |
!** usage: ** |
| 427 |
!** ** |
| 428 |
!** save is called whenever the new verlet list is constructed. ** |
| 429 |
!******************************************************************* |
| 430 |
|
| 431 |
integer :: i |
| 432 |
|
| 433 |
integer :: natoms |
| 434 |
|
| 435 |
double precision, dimension(natoms) :: rx, ry, rz |
| 436 |
double precision, dimension(natoms) :: rx0, ry0, rz0 |
| 437 |
|
| 438 |
!******************************************************************* |
| 439 |
|
| 440 |
do i = 1, natoms |
| 441 |
|
| 442 |
rx0(i) = rx(i) |
| 443 |
ry0(i) = ry(i) |
| 444 |
rz0(i) = rz(i) |
| 445 |
|
| 446 |
end do |
| 447 |
|
| 448 |
return |
| 449 |
end subroutine long_range_save |
| 450 |
|
| 451 |
end module long_range_forces |