| 1 |
module shapes |
| 2 |
|
| 3 |
use force_globals |
| 4 |
use definitions |
| 5 |
use atype_module |
| 6 |
use vector_class |
| 7 |
use simulation |
| 8 |
use status |
| 9 |
#ifdef IS_MPI |
| 10 |
use mpiSimulation |
| 11 |
#endif |
| 12 |
implicit none |
| 13 |
|
| 14 |
PRIVATE |
| 15 |
|
| 16 |
INTEGER, PARAMETER:: CHEBYSHEV_TN = 1 |
| 17 |
INTEGER, PARAMETER:: CHEBYSHEV_UN = 2 |
| 18 |
INTEGER, PARAMETER:: LAGUERRE = 3 |
| 19 |
INTEGER, PARAMETER:: HERMITE = 4 |
| 20 |
INTEGER, PARAMETER:: SH_COS = 0 |
| 21 |
INTEGER, PARAMETER:: SH_SIN = 1 |
| 22 |
|
| 23 |
logical, save :: haveShapeMap = .false. |
| 24 |
|
| 25 |
public :: do_shape_pair |
| 26 |
|
| 27 |
type :: ShapeList |
| 28 |
integer :: nContactFuncs = 0 |
| 29 |
integer :: nRangeFuncs = 0 |
| 30 |
integer :: nStrengthFuncs = 0 |
| 31 |
integer :: bigL = 0 |
| 32 |
integer :: bigM = 0 |
| 33 |
integer, allocatable, dimension(:) :: ContactFuncLValue |
| 34 |
integer, allocatable, dimension(:) :: ContactFuncMValue |
| 35 |
integer, allocatable, dimension(:) :: ContactFunctionType |
| 36 |
real(kind=dp), allocatable, dimension(:) :: ContactFuncCoefficient |
| 37 |
integer, allocatable, dimension(:) :: RangeFuncLValue |
| 38 |
integer, allocatable, dimension(:) :: RangeFuncMValue |
| 39 |
integer, allocatable, dimension(:) :: RangeFunctionType |
| 40 |
real(kind=dp), allocatable, dimension(:) :: RangeFuncCoefficient |
| 41 |
integer, allocatable, dimension(:) :: StrengthFuncLValue |
| 42 |
integer, allocatable, dimension(:) :: StrengthFuncMValue |
| 43 |
integer, allocatable, dimension(:) :: StrengthFunctionType |
| 44 |
real(kind=dp), allocatable, dimension(:) :: StrengthFuncCoefficient |
| 45 |
logical :: isLJ = .false. |
| 46 |
real ( kind = dp ) :: epsilon = 0.0_dp |
| 47 |
real ( kind = dp ) :: sigma = 0.0_dp |
| 48 |
end type ShapeList |
| 49 |
|
| 50 |
type(ShapeList), dimension(:),allocatable :: ShapeMap |
| 51 |
|
| 52 |
integer :: lmax |
| 53 |
real (kind=dp), allocatable, dimension(:,:) :: plm_i, dlm_i, plm_j, dlm_j |
| 54 |
real (kind=dp), allocatable, dimension(:) :: tm_i, dtm_i, um_i, dum_i |
| 55 |
real (kind=dp), allocatable, dimension(:) :: tm_j, dtm_j, um_j, dum_j |
| 56 |
|
| 57 |
contains |
| 58 |
|
| 59 |
subroutine createShapeMap(status) |
| 60 |
integer :: nAtypes |
| 61 |
integer :: status |
| 62 |
integer :: i |
| 63 |
real (kind=DP) :: thisDP |
| 64 |
logical :: thisProperty |
| 65 |
|
| 66 |
status = 0 |
| 67 |
|
| 68 |
nAtypes = getSize(atypes) |
| 69 |
|
| 70 |
if (nAtypes == 0) then |
| 71 |
status = -1 |
| 72 |
return |
| 73 |
end if |
| 74 |
|
| 75 |
if (.not. allocated(ShapeMap)) then |
| 76 |
allocate(ShapeMap(nAtypes)) |
| 77 |
endif |
| 78 |
|
| 79 |
do i = 1, nAtypes |
| 80 |
|
| 81 |
call getElementProperty(atypes, i, "is_SHAPE", thisProperty) |
| 82 |
|
| 83 |
if (thisProperty) then |
| 84 |
|
| 85 |
! do all of the shape stuff |
| 86 |
|
| 87 |
endif |
| 88 |
|
| 89 |
call getElementProperty(atypes, i, "is_LJ", thisProperty) |
| 90 |
|
| 91 |
if (thisProperty) then |
| 92 |
ShapeMap(i)%isLJ = .true. |
| 93 |
call getElementProperty(atypes, i, "lj_epsilon", thisDP) |
| 94 |
ShapeMap(i)%epsilon = thisDP |
| 95 |
call getElementProperty(atypes, i, "lj_sigma", thisDP) |
| 96 |
ShapeMap(i)%sigma = thisDP |
| 97 |
else |
| 98 |
ShapeMap(i)%isLJ = .false. |
| 99 |
endif |
| 100 |
|
| 101 |
|
| 102 |
end do |
| 103 |
|
| 104 |
haveShapeMap = .true. |
| 105 |
|
| 106 |
end subroutine createShapeMap |
| 107 |
|
| 108 |
|
| 109 |
|
| 110 |
subroutine do_shape_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
| 111 |
pot, A, f, t, do_pot) |
| 112 |
|
| 113 |
integer, intent(in) :: atom1, atom2 |
| 114 |
real (kind=dp), intent(inout) :: rij, r2 |
| 115 |
real (kind=dp), dimension(3), intent(in) :: d |
| 116 |
real (kind=dp), dimension(3), intent(inout) :: fpair |
| 117 |
real (kind=dp) :: pot, vpair, sw |
| 118 |
real (kind=dp), dimension(9,nLocal) :: A |
| 119 |
real (kind=dp), dimension(3,nLocal) :: f |
| 120 |
real (kind=dp), dimension(3,nLocal) :: t |
| 121 |
logical, intent(in) :: do_pot |
| 122 |
|
| 123 |
real (kind=dp) :: r3, r5, rt2, rt3, rt5, rt6, rt11, rt12, rt126 |
| 124 |
integer :: me1, me2, l, m, lm, id1, id2, localError, function_type |
| 125 |
real (kind=dp) :: sigma_i, s_i, eps_i, sigma_j, s_j, eps_j |
| 126 |
real (kind=dp) :: coeff |
| 127 |
|
| 128 |
real (kind=dp) :: dsigmaidx, dsigmaidy, dsigmaidz |
| 129 |
real (kind=dp) :: dsigmaidux, dsigmaiduy, dsigmaiduz |
| 130 |
real (kind=dp) :: dsigmajdx, dsigmajdy, dsigmajdz |
| 131 |
real (kind=dp) :: dsigmajdux, dsigmajduy, dsigmajduz |
| 132 |
|
| 133 |
real (kind=dp) :: dsidx, dsidy, dsidz |
| 134 |
real (kind=dp) :: dsidux, dsiduy, dsiduz |
| 135 |
real (kind=dp) :: dsjdx, dsjdy, dsjdz |
| 136 |
real (kind=dp) :: dsjdux, dsjduy, dsjduz |
| 137 |
|
| 138 |
real (kind=dp) :: depsidx, depsidy, depsidz |
| 139 |
real (kind=dp) :: depsidux, depsiduy, depsiduz |
| 140 |
real (kind=dp) :: depsjdx, depsjdy, depsjdz |
| 141 |
real (kind=dp) :: depsjdux, depsjduy, depsjduz |
| 142 |
|
| 143 |
real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2 |
| 144 |
|
| 145 |
real (kind=dp) :: proji, proji3, projj, projj3 |
| 146 |
real (kind=dp) :: cti, ctj, cpi, cpj, spi, spj |
| 147 |
real (kind=dp) :: Phunc, sigma, s, eps, rtdenom, rt |
| 148 |
|
| 149 |
real (kind=dp) :: dctidx, dctidy, dctidz |
| 150 |
real (kind=dp) :: dctidux, dctiduy, dctiduz |
| 151 |
real (kind=dp) :: dctjdx, dctjdy, dctjdz |
| 152 |
real (kind=dp) :: dctjdux, dctjduy, dctjduz |
| 153 |
|
| 154 |
real (kind=dp) :: dcpidx, dcpidy, dcpidz |
| 155 |
real (kind=dp) :: dcpidux, dcpiduy, dcpiduz |
| 156 |
real (kind=dp) :: dcpjdx, dcpjdy, dcpjdz |
| 157 |
real (kind=dp) :: dcpjdux, dcpjduy, dcpjduz |
| 158 |
|
| 159 |
real (kind=dp) :: dspidx, dspidy, dspidz |
| 160 |
real (kind=dp) :: dspidux, dspiduy, dspiduz |
| 161 |
real (kind=dp) :: dspjdx, dspjdy, dspjdz |
| 162 |
real (kind=dp) :: dspjdux, dspjduy, dspjduz |
| 163 |
|
| 164 |
real (kind=dp) :: dPhuncdX, dPhuncdY, dPhuncdZ |
| 165 |
real (kind=dp) :: dPhuncdUx, dPhuncdUy, dPhuncdUz |
| 166 |
|
| 167 |
real (kind=dp) :: dsigmadxi, dsigmadyi, dsigmadzi |
| 168 |
real (kind=dp) :: dsigmaduxi, dsigmaduyi, dsigmaduzi |
| 169 |
real (kind=dp) :: dsigmadxj, dsigmadyj, dsigmadzj |
| 170 |
real (kind=dp) :: dsigmaduxj, dsigmaduyj, dsigmaduzj |
| 171 |
|
| 172 |
real (kind=dp) :: dsdxi, dsdyi, dsdzi |
| 173 |
real (kind=dp) :: dsduxi, dsduyi, dsduzi |
| 174 |
real (kind=dp) :: dsdxj, dsdyj, dsdzj |
| 175 |
real (kind=dp) :: dsduxj, dsduyj, dsduzj |
| 176 |
|
| 177 |
real (kind=dp) :: depsdxi, depsdyi, depsdzi |
| 178 |
real (kind=dp) :: depsduxi, depsduyi, depsduzi |
| 179 |
real (kind=dp) :: depsdxj, depsdyj, depsdzj |
| 180 |
real (kind=dp) :: depsduxj, depsduyj, depsduzj |
| 181 |
|
| 182 |
real (kind=dp) :: drtdxi, drtdyi, drtdzi |
| 183 |
real (kind=dp) :: drtduxi, drtduyi, drtduzi |
| 184 |
real (kind=dp) :: drtdxj, drtdyj, drtdzj |
| 185 |
real (kind=dp) :: drtduxj, drtduyj, drtduzj |
| 186 |
|
| 187 |
real (kind=dp) :: drdxi, drdyi, drdzi |
| 188 |
real (kind=dp) :: drduxi, drduyi, drduzi |
| 189 |
real (kind=dp) :: drdxj, drdyj, drdzj |
| 190 |
real (kind=dp) :: drduxj, drduyj, drduzj |
| 191 |
|
| 192 |
real (kind=dp) :: dvdxi, dvdyi, dvdzi |
| 193 |
real (kind=dp) :: dvduxi, dvduyi, dvduzi |
| 194 |
real (kind=dp) :: dvdxj, dvdyj, dvdzj |
| 195 |
real (kind=dp) :: dvduxj, dvduyj, dvduzj |
| 196 |
|
| 197 |
real (kind=dp) :: fxi, fyi, fzi, fxj, fyj, fzj |
| 198 |
real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj |
| 199 |
real (kind=dp) :: fxii, fyii, fzii, fxij, fyij, fzij |
| 200 |
real (kind=dp) :: fxji, fyji, fzji, fxjj, fyjj, fzjj |
| 201 |
real (kind=dp) :: fxradial, fyradial, fzradial |
| 202 |
|
| 203 |
if (.not.haveShapeMap) then |
| 204 |
localError = 0 |
| 205 |
call createShapeMap(localError) |
| 206 |
if ( localError .ne. 0 ) then |
| 207 |
call handleError("calc_shape", "ShapeMap creation failed!") |
| 208 |
return |
| 209 |
end if |
| 210 |
endif |
| 211 |
|
| 212 |
!! We assume that the rotation matrices have already been calculated |
| 213 |
!! and placed in the A array. |
| 214 |
|
| 215 |
r3 = r2*rij |
| 216 |
r5 = r3*r2 |
| 217 |
|
| 218 |
drdxi = -d(1) / rij |
| 219 |
drdyi = -d(2) / rij |
| 220 |
drdzi = -d(3) / rij |
| 221 |
|
| 222 |
drdxj = d(1) / rij |
| 223 |
drdyj = d(2) / rij |
| 224 |
drdzj = d(3) / rij |
| 225 |
|
| 226 |
#ifdef IS_MPI |
| 227 |
me1 = atid_Row(atom1) |
| 228 |
me2 = atid_Col(atom2) |
| 229 |
#else |
| 230 |
me1 = atid(atom1) |
| 231 |
me2 = atid(atom2) |
| 232 |
#endif |
| 233 |
|
| 234 |
if (ShapeMap(me1)%isLJ) then |
| 235 |
sigma_i = ShapeMap(me1)%sigma |
| 236 |
s_i = ShapeMap(me1)%sigma |
| 237 |
eps_i = ShapeMap(me1)%epsilon |
| 238 |
dsigmaidx = 0.0d0 |
| 239 |
dsigmaidy = 0.0d0 |
| 240 |
dsigmaidz = 0.0d0 |
| 241 |
dsigmaidux = 0.0d0 |
| 242 |
dsigmaiduy = 0.0d0 |
| 243 |
dsigmaiduz = 0.0d0 |
| 244 |
dsidx = 0.0d0 |
| 245 |
dsidy = 0.0d0 |
| 246 |
dsidz = 0.0d0 |
| 247 |
dsidux = 0.0d0 |
| 248 |
dsiduy = 0.0d0 |
| 249 |
dsiduz = 0.0d0 |
| 250 |
depsidx = 0.0d0 |
| 251 |
depsidy = 0.0d0 |
| 252 |
depsidz = 0.0d0 |
| 253 |
depsidux = 0.0d0 |
| 254 |
depsiduy = 0.0d0 |
| 255 |
depsiduz = 0.0d0 |
| 256 |
else |
| 257 |
|
| 258 |
#ifdef IS_MPI |
| 259 |
! rotate the inter-particle separation into the two different |
| 260 |
! body-fixed coordinate systems: |
| 261 |
|
| 262 |
xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3) |
| 263 |
yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3) |
| 264 |
zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3) |
| 265 |
|
| 266 |
#else |
| 267 |
! rotate the inter-particle separation into the two different |
| 268 |
! body-fixed coordinate systems: |
| 269 |
|
| 270 |
xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3) |
| 271 |
yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3) |
| 272 |
zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3) |
| 273 |
|
| 274 |
#endif |
| 275 |
|
| 276 |
xi2 = xi*xi |
| 277 |
yi2 = yi*yi |
| 278 |
zi2 = zi*zi |
| 279 |
|
| 280 |
proji = sqrt(xi2 + yi2) |
| 281 |
proji3 = proji*proji*proji |
| 282 |
|
| 283 |
cti = zi / rij |
| 284 |
dctidx = - zi * xi / r3 |
| 285 |
dctidy = - zi * yi / r3 |
| 286 |
dctidz = 1.0d0 / rij - zi2 / r3 |
| 287 |
dctidux = yi / rij |
| 288 |
dctiduy = -xi / rij |
| 289 |
dctiduz = 0.0d0 |
| 290 |
|
| 291 |
cpi = xi / proji |
| 292 |
dcpidx = 1.0d0 / proji - xi2 / proji3 |
| 293 |
dcpidy = - xi * yi / proji3 |
| 294 |
dcpidz = 0.0d0 |
| 295 |
dcpidux = xi * yi * zi / proji3 |
| 296 |
dcpiduy = -zi * (1.0d0 / proji - xi2 / proji3) |
| 297 |
dcpiduz = -yi * (1.0d0 / proji - xi2 / proji3) - (xi2 * yi / proji3) |
| 298 |
|
| 299 |
spi = yi / proji |
| 300 |
dspidx = - xi * yi / proji3 |
| 301 |
dspidy = 1.0d0 / proji - yi2 / proji3 |
| 302 |
dspidz = 0.0d0 |
| 303 |
dspidux = -zi * (1.0d0 / proji - yi2 / proji3) |
| 304 |
dspiduy = xi * yi * zi / proji3 |
| 305 |
dspiduz = xi * (1.0d0 / proji - yi2 / proji3) + (xi * yi2 / proji3) |
| 306 |
|
| 307 |
call Associated_Legendre(cti, ShapeMap(me1)%bigL, & |
| 308 |
ShapeMap(me1)%bigM, lmax, plm_i, dlm_i) |
| 309 |
|
| 310 |
call Orthogonal_Polynomial(cpi, ShapeMap(me1)%bigM, CHEBYSHEV_TN, & |
| 311 |
tm_i, dtm_i) |
| 312 |
call Orthogonal_Polynomial(cpi, ShapeMap(me1)%bigM, CHEBYSHEV_UN, & |
| 313 |
um_i, dum_i) |
| 314 |
|
| 315 |
sigma_i = 0.0d0 |
| 316 |
s_i = 0.0d0 |
| 317 |
eps_i = 0.0d0 |
| 318 |
dsigmaidx = 0.0d0 |
| 319 |
dsigmaidy = 0.0d0 |
| 320 |
dsigmaidz = 0.0d0 |
| 321 |
dsigmaidux = 0.0d0 |
| 322 |
dsigmaiduy = 0.0d0 |
| 323 |
dsigmaiduz = 0.0d0 |
| 324 |
dsidx = 0.0d0 |
| 325 |
dsidy = 0.0d0 |
| 326 |
dsidz = 0.0d0 |
| 327 |
dsidux = 0.0d0 |
| 328 |
dsiduy = 0.0d0 |
| 329 |
dsiduz = 0.0d0 |
| 330 |
depsidx = 0.0d0 |
| 331 |
depsidy = 0.0d0 |
| 332 |
depsidz = 0.0d0 |
| 333 |
depsidux = 0.0d0 |
| 334 |
depsiduy = 0.0d0 |
| 335 |
depsiduz = 0.0d0 |
| 336 |
|
| 337 |
do lm = 1, ShapeMap(me1)%nContactFuncs |
| 338 |
l = ShapeMap(me1)%ContactFuncLValue(lm) |
| 339 |
m = ShapeMap(me1)%ContactFuncMValue(lm) |
| 340 |
coeff = ShapeMap(me1)%ContactFuncCoefficient(lm) |
| 341 |
function_type = ShapeMap(me1)%ContactFunctionType(lm) |
| 342 |
|
| 343 |
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 344 |
Phunc = coeff * tm_i(m) |
| 345 |
dPhuncdX = coeff * dtm_i(m) * dcpidx |
| 346 |
dPhuncdY = coeff * dtm_i(m) * dcpidy |
| 347 |
dPhuncdZ = coeff * dtm_i(m) * dcpidz |
| 348 |
dPhuncdUz = coeff * dtm_i(m) * dcpidux |
| 349 |
dPhuncdUy = coeff * dtm_i(m) * dcpiduy |
| 350 |
dPhuncdUz = coeff * dtm_i(m) * dcpiduz |
| 351 |
else |
| 352 |
Phunc = coeff * spi * um_i(m-1) |
| 353 |
dPhuncdX = coeff * (spi * dum_i(m-1) * dcpidx + dspidx *um_i(m-1)) |
| 354 |
dPhuncdY = coeff * (spi * dum_i(m-1) * dcpidy + dspidy *um_i(m-1)) |
| 355 |
dPhuncdZ = coeff * (spi * dum_i(m-1) * dcpidz + dspidz *um_i(m-1)) |
| 356 |
dPhuncdUx = coeff*(spi * dum_i(m-1)*dcpidux + dspidux *um_i(m-1)) |
| 357 |
dPhuncdUy = coeff*(spi * dum_i(m-1)*dcpiduy + dspiduy *um_i(m-1)) |
| 358 |
dPhuncdUz = coeff*(spi * dum_i(m-1)*dcpiduz + dspiduz *um_i(m-1)) |
| 359 |
endif |
| 360 |
|
| 361 |
sigma_i = sigma_i + plm_i(l,m)*Phunc |
| 362 |
|
| 363 |
dsigmaidx = dsigmaidx + plm_i(l,m)*dPhuncdX + & |
| 364 |
Phunc * dlm_i(l,m) * dctidx |
| 365 |
dsigmaidy = dsigmaidy + plm_i(l,m)*dPhuncdY + & |
| 366 |
Phunc * dlm_i(l,m) * dctidy |
| 367 |
dsigmaidz = dsigmaidz + plm_i(l,m)*dPhuncdZ + & |
| 368 |
Phunc * dlm_i(l,m) * dctidz |
| 369 |
|
| 370 |
dsigmaidux = dsigmaidux + plm_i(l,m)* dPhuncdUx + & |
| 371 |
Phunc * dlm_i(l,m) * dctidux |
| 372 |
dsigmaiduy = dsigmaiduy + plm_i(l,m)* dPhuncdUy + & |
| 373 |
Phunc * dlm_i(l,m) * dctiduy |
| 374 |
dsigmaiduz = dsigmaiduz + plm_i(l,m)* dPhuncdUz + & |
| 375 |
Phunc * dlm_i(l,m) * dctiduz |
| 376 |
|
| 377 |
end do |
| 378 |
|
| 379 |
do lm = 1, ShapeMap(me1)%nRangeFuncs |
| 380 |
l = ShapeMap(me1)%RangeFuncLValue(lm) |
| 381 |
m = ShapeMap(me1)%RangeFuncMValue(lm) |
| 382 |
coeff = ShapeMap(me1)%RangeFuncCoefficient(lm) |
| 383 |
function_type = ShapeMap(me1)%RangeFunctionType(lm) |
| 384 |
|
| 385 |
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 386 |
Phunc = coeff * tm_i(m) |
| 387 |
dPhuncdX = coeff * dtm_i(m) * dcpidx |
| 388 |
dPhuncdY = coeff * dtm_i(m) * dcpidy |
| 389 |
dPhuncdZ = coeff * dtm_i(m) * dcpidz |
| 390 |
dPhuncdUz = coeff * dtm_i(m) * dcpidux |
| 391 |
dPhuncdUy = coeff * dtm_i(m) * dcpiduy |
| 392 |
dPhuncdUz = coeff * dtm_i(m) * dcpiduz |
| 393 |
else |
| 394 |
Phunc = coeff * spi * um_i(m-1) |
| 395 |
dPhuncdX = coeff * (spi * dum_i(m-1) * dcpidx + dspidx *um_i(m-1)) |
| 396 |
dPhuncdY = coeff * (spi * dum_i(m-1) * dcpidy + dspidy *um_i(m-1)) |
| 397 |
dPhuncdZ = coeff * (spi * dum_i(m-1) * dcpidz + dspidz *um_i(m-1)) |
| 398 |
dPhuncdUx = coeff*(spi * dum_i(m-1)*dcpidux + dspidux *um_i(m-1)) |
| 399 |
dPhuncdUy = coeff*(spi * dum_i(m-1)*dcpiduy + dspiduy *um_i(m-1)) |
| 400 |
dPhuncdUz = coeff*(spi * dum_i(m-1)*dcpiduz + dspiduz *um_i(m-1)) |
| 401 |
endif |
| 402 |
|
| 403 |
s_i = s_i + plm_i(l,m)*Phunc |
| 404 |
|
| 405 |
dsidx = dsidx + plm_i(l,m)*dPhuncdX + & |
| 406 |
Phunc * dlm_i(l,m) * dctidx |
| 407 |
dsidy = dsidy + plm_i(l,m)*dPhuncdY + & |
| 408 |
Phunc * dlm_i(l,m) * dctidy |
| 409 |
dsidz = dsidz + plm_i(l,m)*dPhuncdZ + & |
| 410 |
Phunc * dlm_i(l,m) * dctidz |
| 411 |
|
| 412 |
dsidux = dsidux + plm_i(l,m)* dPhuncdUx + & |
| 413 |
Phunc * dlm_i(l,m) * dctidux |
| 414 |
dsiduy = dsiduy + plm_i(l,m)* dPhuncdUy + & |
| 415 |
Phunc * dlm_i(l,m) * dctiduy |
| 416 |
dsiduz = dsiduz + plm_i(l,m)* dPhuncdUz + & |
| 417 |
Phunc * dlm_i(l,m) * dctiduz |
| 418 |
|
| 419 |
end do |
| 420 |
|
| 421 |
do lm = 1, ShapeMap(me1)%nStrengthFuncs |
| 422 |
l = ShapeMap(me1)%StrengthFuncLValue(lm) |
| 423 |
m = ShapeMap(me1)%StrengthFuncMValue(lm) |
| 424 |
coeff = ShapeMap(me1)%StrengthFuncCoefficient(lm) |
| 425 |
function_type = ShapeMap(me1)%StrengthFunctionType(lm) |
| 426 |
|
| 427 |
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 428 |
Phunc = coeff * tm_i(m) |
| 429 |
dPhuncdX = coeff * dtm_i(m) * dcpidx |
| 430 |
dPhuncdY = coeff * dtm_i(m) * dcpidy |
| 431 |
dPhuncdZ = coeff * dtm_i(m) * dcpidz |
| 432 |
dPhuncdUz = coeff * dtm_i(m) * dcpidux |
| 433 |
dPhuncdUy = coeff * dtm_i(m) * dcpiduy |
| 434 |
dPhuncdUz = coeff * dtm_i(m) * dcpiduz |
| 435 |
else |
| 436 |
Phunc = coeff * spi * um_i(m-1) |
| 437 |
dPhuncdX = coeff * (spi * dum_i(m-1) * dcpidx + dspidx *um_i(m-1)) |
| 438 |
dPhuncdY = coeff * (spi * dum_i(m-1) * dcpidy + dspidy *um_i(m-1)) |
| 439 |
dPhuncdZ = coeff * (spi * dum_i(m-1) * dcpidz + dspidz *um_i(m-1)) |
| 440 |
dPhuncdUx = coeff*(spi * dum_i(m-1)*dcpidux + dspidux *um_i(m-1)) |
| 441 |
dPhuncdUy = coeff*(spi * dum_i(m-1)*dcpiduy + dspiduy *um_i(m-1)) |
| 442 |
dPhuncdUz = coeff*(spi * dum_i(m-1)*dcpiduz + dspiduz *um_i(m-1)) |
| 443 |
endif |
| 444 |
|
| 445 |
eps_i = eps_i + plm_i(l,m)*Phunc |
| 446 |
|
| 447 |
depsidx = depsidx + plm_i(l,m)*dPhuncdX + & |
| 448 |
Phunc * dlm_i(l,m) * dctidx |
| 449 |
depsidy = depsidy + plm_i(l,m)*dPhuncdY + & |
| 450 |
Phunc * dlm_i(l,m) * dctidy |
| 451 |
depsidz = depsidz + plm_i(l,m)*dPhuncdZ + & |
| 452 |
Phunc * dlm_i(l,m) * dctidz |
| 453 |
|
| 454 |
depsidux = depsidux + plm_i(l,m)* dPhuncdUx + & |
| 455 |
Phunc * dlm_i(l,m) * dctidux |
| 456 |
depsiduy = depsiduy + plm_i(l,m)* dPhuncdUy + & |
| 457 |
Phunc * dlm_i(l,m) * dctiduy |
| 458 |
depsiduz = depsiduz + plm_i(l,m)* dPhuncdUz + & |
| 459 |
Phunc * dlm_i(l,m) * dctiduz |
| 460 |
|
| 461 |
end do |
| 462 |
|
| 463 |
endif |
| 464 |
|
| 465 |
! now do j: |
| 466 |
|
| 467 |
if (ShapeMap(me2)%isLJ) then |
| 468 |
sigma_j = ShapeMap(me2)%sigma |
| 469 |
s_j = ShapeMap(me2)%sigma |
| 470 |
eps_j = ShapeMap(me2)%epsilon |
| 471 |
dsigmajdx = 0.0d0 |
| 472 |
dsigmajdy = 0.0d0 |
| 473 |
dsigmajdz = 0.0d0 |
| 474 |
dsigmajdux = 0.0d0 |
| 475 |
dsigmajduy = 0.0d0 |
| 476 |
dsigmajduz = 0.0d0 |
| 477 |
dsjdx = 0.0d0 |
| 478 |
dsjdy = 0.0d0 |
| 479 |
dsjdz = 0.0d0 |
| 480 |
dsjdux = 0.0d0 |
| 481 |
dsjduy = 0.0d0 |
| 482 |
dsjduz = 0.0d0 |
| 483 |
depsjdx = 0.0d0 |
| 484 |
depsjdy = 0.0d0 |
| 485 |
depsjdz = 0.0d0 |
| 486 |
depsjdux = 0.0d0 |
| 487 |
depsjduy = 0.0d0 |
| 488 |
depsjduz = 0.0d0 |
| 489 |
else |
| 490 |
|
| 491 |
#ifdef IS_MPI |
| 492 |
! rotate the inter-particle separation into the two different |
| 493 |
! body-fixed coordinate systems: |
| 494 |
! negative sign because this is the vector from j to i: |
| 495 |
|
| 496 |
xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3)) |
| 497 |
yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3)) |
| 498 |
zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3)) |
| 499 |
#else |
| 500 |
! rotate the inter-particle separation into the two different |
| 501 |
! body-fixed coordinate systems: |
| 502 |
! negative sign because this is the vector from j to i: |
| 503 |
|
| 504 |
xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3)) |
| 505 |
yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3)) |
| 506 |
zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3)) |
| 507 |
#endif |
| 508 |
|
| 509 |
xj2 = xj*xj |
| 510 |
yj2 = yj*yj |
| 511 |
zj2 = zj*zj |
| 512 |
|
| 513 |
projj = sqrt(xj2 + yj2) |
| 514 |
projj3 = projj*projj*projj |
| 515 |
|
| 516 |
ctj = zj / rij |
| 517 |
dctjdx = - zj * xj / r3 |
| 518 |
dctjdy = - zj * yj / r3 |
| 519 |
dctjdz = 1.0d0 / rij - zj2 / r3 |
| 520 |
dctjdux = yj / rij |
| 521 |
dctjduy = -xj / rij |
| 522 |
dctjduz = 0.0d0 |
| 523 |
|
| 524 |
cpj = xj / projj |
| 525 |
dcpjdx = 1.0d0 / projj - xj2 / projj3 |
| 526 |
dcpjdy = - xj * yj / projj3 |
| 527 |
dcpjdz = 0.0d0 |
| 528 |
dcpjdux = xj * yj * zj / projj3 |
| 529 |
dcpjduy = -zj * (1.0d0 / projj - xj2 / projj3) |
| 530 |
dcpjduz = -yj * (1.0d0 / projj - xj2 / projj3) - (xj2 * yj / projj3) |
| 531 |
|
| 532 |
spj = yj / projj |
| 533 |
dspjdx = - xj * yj / projj3 |
| 534 |
dspjdy = 1.0d0 / projj - yj2 / projj3 |
| 535 |
dspjdz = 0.0d0 |
| 536 |
dspjdux = -zj * (1.0d0 / projj - yj2 / projj3) |
| 537 |
dspjduy = xj * yj * zj / projj3 |
| 538 |
dspjduz = xj * (1.0d0 / projj - yi2 / projj3) + (xj * yj2 / projj3) |
| 539 |
|
| 540 |
call Associated_Legendre(ctj, ShapeMap(me2)%bigL, & |
| 541 |
ShapeMap(me2)%bigM, lmax, plm_j, dlm_j) |
| 542 |
|
| 543 |
call Orthogonal_Polynomial(cpj, ShapeMap(me2)%bigM, CHEBYSHEV_TN, & |
| 544 |
tm_j, dtm_j) |
| 545 |
call Orthogonal_Polynomial(cpj, ShapeMap(me2)%bigM, CHEBYSHEV_UN, & |
| 546 |
um_j, dum_j) |
| 547 |
|
| 548 |
sigma_j = 0.0d0 |
| 549 |
s_j = 0.0d0 |
| 550 |
eps_j = 0.0d0 |
| 551 |
dsigmajdx = 0.0d0 |
| 552 |
dsigmajdy = 0.0d0 |
| 553 |
dsigmajdz = 0.0d0 |
| 554 |
dsigmajdux = 0.0d0 |
| 555 |
dsigmajduy = 0.0d0 |
| 556 |
dsigmajduz = 0.0d0 |
| 557 |
dsjdx = 0.0d0 |
| 558 |
dsjdy = 0.0d0 |
| 559 |
dsjdz = 0.0d0 |
| 560 |
dsjdux = 0.0d0 |
| 561 |
dsjduy = 0.0d0 |
| 562 |
dsjduz = 0.0d0 |
| 563 |
depsjdx = 0.0d0 |
| 564 |
depsjdy = 0.0d0 |
| 565 |
depsjdz = 0.0d0 |
| 566 |
depsjdux = 0.0d0 |
| 567 |
depsjduy = 0.0d0 |
| 568 |
depsjduz = 0.0d0 |
| 569 |
|
| 570 |
do lm = 1, ShapeMap(me2)%nContactFuncs |
| 571 |
l = ShapeMap(me2)%ContactFuncLValue(lm) |
| 572 |
m = ShapeMap(me2)%ContactFuncMValue(lm) |
| 573 |
coeff = ShapeMap(me2)%ContactFuncCoefficient(lm) |
| 574 |
function_type = ShapeMap(me2)%ContactFunctionType(lm) |
| 575 |
|
| 576 |
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 577 |
Phunc = coeff * tm_j(m) |
| 578 |
dPhuncdX = coeff * dtm_j(m) * dcpjdx |
| 579 |
dPhuncdY = coeff * dtm_j(m) * dcpjdy |
| 580 |
dPhuncdZ = coeff * dtm_j(m) * dcpjdz |
| 581 |
dPhuncdUz = coeff * dtm_j(m) * dcpjdux |
| 582 |
dPhuncdUy = coeff * dtm_j(m) * dcpjduy |
| 583 |
dPhuncdUz = coeff * dtm_j(m) * dcpjduz |
| 584 |
else |
| 585 |
Phunc = coeff * spj * um_j(m-1) |
| 586 |
dPhuncdX = coeff * (spj * dum_j(m-1) * dcpjdx + dspjdx *um_j(m-1)) |
| 587 |
dPhuncdY = coeff * (spj * dum_j(m-1) * dcpjdy + dspjdy *um_j(m-1)) |
| 588 |
dPhuncdZ = coeff * (spj * dum_j(m-1) * dcpjdz + dspjdz *um_j(m-1)) |
| 589 |
dPhuncdUx = coeff*(spj * dum_j(m-1)*dcpjdux + dspjdux *um_j(m-1)) |
| 590 |
dPhuncdUy = coeff*(spj * dum_j(m-1)*dcpjduy + dspjduy *um_j(m-1)) |
| 591 |
dPhuncdUz = coeff*(spj * dum_j(m-1)*dcpjduz + dspjduz *um_j(m-1)) |
| 592 |
endif |
| 593 |
|
| 594 |
sigma_j = sigma_j + plm_j(l,m)*Phunc |
| 595 |
|
| 596 |
dsigmajdx = dsigmajdx + plm_j(l,m)*dPhuncdX + & |
| 597 |
Phunc * dlm_j(l,m) * dctjdx |
| 598 |
dsigmajdy = dsigmajdy + plm_j(l,m)*dPhuncdY + & |
| 599 |
Phunc * dlm_j(l,m) * dctjdy |
| 600 |
dsigmajdz = dsigmajdz + plm_j(l,m)*dPhuncdZ + & |
| 601 |
Phunc * dlm_j(l,m) * dctjdz |
| 602 |
|
| 603 |
dsigmajdux = dsigmajdux + plm_j(l,m)* dPhuncdUx + & |
| 604 |
Phunc * dlm_j(l,m) * dctjdux |
| 605 |
dsigmajduy = dsigmajduy + plm_j(l,m)* dPhuncdUy + & |
| 606 |
Phunc * dlm_j(l,m) * dctjduy |
| 607 |
dsigmajduz = dsigmajduz + plm_j(l,m)* dPhuncdUz + & |
| 608 |
Phunc * dlm_j(l,m) * dctjduz |
| 609 |
|
| 610 |
end do |
| 611 |
|
| 612 |
do lm = 1, ShapeMap(me2)%nRangeFuncs |
| 613 |
l = ShapeMap(me2)%RangeFuncLValue(lm) |
| 614 |
m = ShapeMap(me2)%RangeFuncMValue(lm) |
| 615 |
coeff = ShapeMap(me2)%RangeFuncCoefficient(lm) |
| 616 |
function_type = ShapeMap(me2)%RangeFunctionType(lm) |
| 617 |
|
| 618 |
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 619 |
Phunc = coeff * tm_j(m) |
| 620 |
dPhuncdX = coeff * dtm_j(m) * dcpjdx |
| 621 |
dPhuncdY = coeff * dtm_j(m) * dcpjdy |
| 622 |
dPhuncdZ = coeff * dtm_j(m) * dcpjdz |
| 623 |
dPhuncdUz = coeff * dtm_j(m) * dcpjdux |
| 624 |
dPhuncdUy = coeff * dtm_j(m) * dcpjduy |
| 625 |
dPhuncdUz = coeff * dtm_j(m) * dcpjduz |
| 626 |
else |
| 627 |
Phunc = coeff * spj * um_j(m-1) |
| 628 |
dPhuncdX = coeff * (spj * dum_j(m-1) * dcpjdx + dspjdx *um_j(m-1)) |
| 629 |
dPhuncdY = coeff * (spj * dum_j(m-1) * dcpjdy + dspjdy *um_j(m-1)) |
| 630 |
dPhuncdZ = coeff * (spj * dum_j(m-1) * dcpjdz + dspjdz *um_j(m-1)) |
| 631 |
dPhuncdUx = coeff*(spj * dum_j(m-1)*dcpjdux + dspjdux *um_j(m-1)) |
| 632 |
dPhuncdUy = coeff*(spj * dum_j(m-1)*dcpjduy + dspjduy *um_j(m-1)) |
| 633 |
dPhuncdUz = coeff*(spj * dum_j(m-1)*dcpjduz + dspjduz *um_j(m-1)) |
| 634 |
endif |
| 635 |
|
| 636 |
s_j = s_j + plm_j(l,m)*Phunc |
| 637 |
|
| 638 |
dsjdx = dsjdx + plm_j(l,m)*dPhuncdX + & |
| 639 |
Phunc * dlm_j(l,m) * dctjdx |
| 640 |
dsjdy = dsjdy + plm_j(l,m)*dPhuncdY + & |
| 641 |
Phunc * dlm_j(l,m) * dctjdy |
| 642 |
dsjdz = dsjdz + plm_j(l,m)*dPhuncdZ + & |
| 643 |
Phunc * dlm_j(l,m) * dctjdz |
| 644 |
|
| 645 |
dsjdux = dsjdux + plm_j(l,m)* dPhuncdUx + & |
| 646 |
Phunc * dlm_j(l,m) * dctjdux |
| 647 |
dsjduy = dsjduy + plm_j(l,m)* dPhuncdUy + & |
| 648 |
Phunc * dlm_j(l,m) * dctjduy |
| 649 |
dsjduz = dsjduz + plm_j(l,m)* dPhuncdUz + & |
| 650 |
Phunc * dlm_j(l,m) * dctjduz |
| 651 |
|
| 652 |
end do |
| 653 |
|
| 654 |
do lm = 1, ShapeMap(me2)%nStrengthFuncs |
| 655 |
l = ShapeMap(me2)%StrengthFuncLValue(lm) |
| 656 |
m = ShapeMap(me2)%StrengthFuncMValue(lm) |
| 657 |
coeff = ShapeMap(me2)%StrengthFuncCoefficient(lm) |
| 658 |
function_type = ShapeMap(me2)%StrengthFunctionType(lm) |
| 659 |
|
| 660 |
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 661 |
Phunc = coeff * tm_j(m) |
| 662 |
dPhuncdX = coeff * dtm_j(m) * dcpjdx |
| 663 |
dPhuncdY = coeff * dtm_j(m) * dcpjdy |
| 664 |
dPhuncdZ = coeff * dtm_j(m) * dcpjdz |
| 665 |
dPhuncdUz = coeff * dtm_j(m) * dcpjdux |
| 666 |
dPhuncdUy = coeff * dtm_j(m) * dcpjduy |
| 667 |
dPhuncdUz = coeff * dtm_j(m) * dcpjduz |
| 668 |
else |
| 669 |
Phunc = coeff * spj * um_j(m-1) |
| 670 |
dPhuncdX = coeff * (spj * dum_j(m-1) * dcpjdx + dspjdx *um_j(m-1)) |
| 671 |
dPhuncdY = coeff * (spj * dum_j(m-1) * dcpjdy + dspjdy *um_j(m-1)) |
| 672 |
dPhuncdZ = coeff * (spj * dum_j(m-1) * dcpjdz + dspjdz *um_j(m-1)) |
| 673 |
dPhuncdUx = coeff*(spj * dum_j(m-1)*dcpjdux + dspjdux *um_j(m-1)) |
| 674 |
dPhuncdUy = coeff*(spj * dum_j(m-1)*dcpjduy + dspjduy *um_j(m-1)) |
| 675 |
dPhuncdUz = coeff*(spj * dum_j(m-1)*dcpjduz + dspjduz *um_j(m-1)) |
| 676 |
endif |
| 677 |
|
| 678 |
eps_j = eps_j + plm_j(l,m)*Phunc |
| 679 |
|
| 680 |
depsjdx = depsjdx + plm_j(l,m)*dPhuncdX + & |
| 681 |
Phunc * dlm_j(l,m) * dctjdx |
| 682 |
depsjdy = depsjdy + plm_j(l,m)*dPhuncdY + & |
| 683 |
Phunc * dlm_j(l,m) * dctjdy |
| 684 |
depsjdz = depsjdz + plm_j(l,m)*dPhuncdZ + & |
| 685 |
Phunc * dlm_j(l,m) * dctjdz |
| 686 |
|
| 687 |
depsjdux = depsjdux + plm_j(l,m)* dPhuncdUx + & |
| 688 |
Phunc * dlm_j(l,m) * dctjdux |
| 689 |
depsjduy = depsjduy + plm_j(l,m)* dPhuncdUy + & |
| 690 |
Phunc * dlm_j(l,m) * dctjduy |
| 691 |
depsjduz = depsjduz + plm_j(l,m)* dPhuncdUz + & |
| 692 |
Phunc * dlm_j(l,m) * dctjduz |
| 693 |
|
| 694 |
end do |
| 695 |
|
| 696 |
endif |
| 697 |
|
| 698 |
! phew, now let's assemble the potential energy: |
| 699 |
|
| 700 |
sigma = 0.5*(sigma_i + sigma_j) |
| 701 |
|
| 702 |
dsigmadxi = 0.5*dsigmaidx |
| 703 |
dsigmadyi = 0.5*dsigmaidy |
| 704 |
dsigmadzi = 0.5*dsigmaidz |
| 705 |
dsigmaduxi = 0.5*dsigmaidux |
| 706 |
dsigmaduyi = 0.5*dsigmaiduy |
| 707 |
dsigmaduzi = 0.5*dsigmaiduz |
| 708 |
|
| 709 |
dsigmadxj = 0.5*dsigmajdx |
| 710 |
dsigmadyj = 0.5*dsigmajdy |
| 711 |
dsigmadzj = 0.5*dsigmajdz |
| 712 |
dsigmaduxj = 0.5*dsigmajdux |
| 713 |
dsigmaduyj = 0.5*dsigmajduy |
| 714 |
dsigmaduzj = 0.5*dsigmajduz |
| 715 |
|
| 716 |
s = 0.5*(s_i + s_j) |
| 717 |
|
| 718 |
dsdxi = 0.5*dsidx |
| 719 |
dsdyi = 0.5*dsidy |
| 720 |
dsdzi = 0.5*dsidz |
| 721 |
dsduxi = 0.5*dsidux |
| 722 |
dsduyi = 0.5*dsiduy |
| 723 |
dsduzi = 0.5*dsiduz |
| 724 |
|
| 725 |
dsdxj = 0.5*dsjdx |
| 726 |
dsdyj = 0.5*dsjdy |
| 727 |
dsdzj = 0.5*dsjdz |
| 728 |
dsduxj = 0.5*dsjdux |
| 729 |
dsduyj = 0.5*dsjduy |
| 730 |
dsduzj = 0.5*dsjduz |
| 731 |
|
| 732 |
eps = sqrt(eps_i * eps_j) |
| 733 |
|
| 734 |
depsdxi = eps_j * depsidx / (2.0d0 * eps) |
| 735 |
depsdyi = eps_j * depsidy / (2.0d0 * eps) |
| 736 |
depsdzi = eps_j * depsidz / (2.0d0 * eps) |
| 737 |
depsduxi = eps_j * depsidux / (2.0d0 * eps) |
| 738 |
depsduyi = eps_j * depsiduy / (2.0d0 * eps) |
| 739 |
depsduzi = eps_j * depsiduz / (2.0d0 * eps) |
| 740 |
|
| 741 |
depsdxj = eps_i * depsjdx / (2.0d0 * eps) |
| 742 |
depsdyj = eps_i * depsjdy / (2.0d0 * eps) |
| 743 |
depsdzj = eps_i * depsjdz / (2.0d0 * eps) |
| 744 |
depsduxj = eps_i * depsjdux / (2.0d0 * eps) |
| 745 |
depsduyj = eps_i * depsjduy / (2.0d0 * eps) |
| 746 |
depsduzj = eps_i * depsjduz / (2.0d0 * eps) |
| 747 |
|
| 748 |
rtdenom = rij-sigma+s |
| 749 |
rt = s / rtdenom |
| 750 |
|
| 751 |
drtdxi = (dsdxi + rt * (drdxi - dsigmadxi + dsdxi)) / rtdenom |
| 752 |
drtdyi = (dsdyi + rt * (drdyi - dsigmadyi + dsdyi)) / rtdenom |
| 753 |
drtdzi = (dsdzi + rt * (drdzi - dsigmadzi + dsdzi)) / rtdenom |
| 754 |
drtduxi = (dsduxi + rt * (drduxi - dsigmaduxi + dsduxi)) / rtdenom |
| 755 |
drtduyi = (dsduyi + rt * (drduyi - dsigmaduyi + dsduyi)) / rtdenom |
| 756 |
drtduzi = (dsduzi + rt * (drduzi - dsigmaduzi + dsduzi)) / rtdenom |
| 757 |
drtdxj = (dsdxj + rt * (drdxj - dsigmadxj + dsdxj)) / rtdenom |
| 758 |
drtdyj = (dsdyj + rt * (drdyj - dsigmadyj + dsdyj)) / rtdenom |
| 759 |
drtdzj = (dsdzj + rt * (drdzj - dsigmadzj + dsdzj)) / rtdenom |
| 760 |
drtduxj = (dsduxj + rt * (drduxj - dsigmaduxj + dsduxj)) / rtdenom |
| 761 |
drtduyj = (dsduyj + rt * (drduyj - dsigmaduyj + dsduyj)) / rtdenom |
| 762 |
drtduzj = (dsduzj + rt * (drduzj - dsigmaduzj + dsduzj)) / rtdenom |
| 763 |
|
| 764 |
rt2 = rt*rt |
| 765 |
rt3 = rt2*rt |
| 766 |
rt5 = rt2*rt3 |
| 767 |
rt6 = rt3*rt3 |
| 768 |
rt11 = rt5*rt6 |
| 769 |
rt12 = rt6*rt6 |
| 770 |
rt126 = rt12 - rt6 |
| 771 |
|
| 772 |
if (do_pot) then |
| 773 |
#ifdef IS_MPI |
| 774 |
pot_row(atom1) = pot_row(atom1) + 2.0d0*eps*rt126*sw |
| 775 |
pot_col(atom2) = pot_col(atom2) + 2.0d0*eps*rt126*sw |
| 776 |
#else |
| 777 |
pot = pot + 4.0d0*eps*rt126*sw |
| 778 |
#endif |
| 779 |
endif |
| 780 |
|
| 781 |
dvdxi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdxi + 4.0d0*depsdxi*rt126 |
| 782 |
dvdyi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdyi + 4.0d0*depsdyi*rt126 |
| 783 |
dvdzi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdzi + 4.0d0*depsdzi*rt126 |
| 784 |
dvduxi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduxi + 4.0d0*depsduxi*rt126 |
| 785 |
dvduyi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduyi + 4.0d0*depsduyi*rt126 |
| 786 |
dvduzi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduzi + 4.0d0*depsduzi*rt126 |
| 787 |
|
| 788 |
dvdxj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdxj + 4.0d0*depsdxj*rt126 |
| 789 |
dvdyj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdyj + 4.0d0*depsdyj*rt126 |
| 790 |
dvdzj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdzj + 4.0d0*depsdzj*rt126 |
| 791 |
dvduxj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduxj + 4.0d0*depsduxj*rt126 |
| 792 |
dvduyj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduyj + 4.0d0*depsduyj*rt126 |
| 793 |
dvduzj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduzj + 4.0d0*depsduzj*rt126 |
| 794 |
|
| 795 |
! do the torques first since they are easy: |
| 796 |
! remember that these are still in the body fixed axes |
| 797 |
|
| 798 |
txi = dvduxi * sw |
| 799 |
tyi = dvduyi * sw |
| 800 |
tzi = dvduzi * sw |
| 801 |
|
| 802 |
txj = dvduxj * sw |
| 803 |
tyj = dvduyj * sw |
| 804 |
tzj = dvduzj * sw |
| 805 |
|
| 806 |
! go back to lab frame using transpose of rotation matrix: |
| 807 |
|
| 808 |
#ifdef IS_MPI |
| 809 |
t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + & |
| 810 |
a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi |
| 811 |
t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + & |
| 812 |
a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi |
| 813 |
t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + & |
| 814 |
a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi |
| 815 |
|
| 816 |
t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + & |
| 817 |
a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj |
| 818 |
t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + & |
| 819 |
a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj |
| 820 |
t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + & |
| 821 |
a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj |
| 822 |
#else |
| 823 |
t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi |
| 824 |
t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi |
| 825 |
t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi |
| 826 |
|
| 827 |
t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj |
| 828 |
t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj |
| 829 |
t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj |
| 830 |
#endif |
| 831 |
! Now, on to the forces: |
| 832 |
|
| 833 |
! first rotate the i terms back into the lab frame: |
| 834 |
|
| 835 |
fxi = dvdxi * sw |
| 836 |
fyi = dvdyi * sw |
| 837 |
fzi = dvdzi * sw |
| 838 |
|
| 839 |
fxj = dvdxj * sw |
| 840 |
fyj = dvdyj * sw |
| 841 |
fzj = dvdzj * sw |
| 842 |
|
| 843 |
#ifdef IS_MPI |
| 844 |
fxii = a_Row(1,atom1)*fxi + a_Row(4,atom1)*fyi + a_Row(7,atom1)*fzi |
| 845 |
fyii = a_Row(2,atom1)*fxi + a_Row(5,atom1)*fyi + a_Row(8,atom1)*fzi |
| 846 |
fzii = a_Row(3,atom1)*fxi + a_Row(6,atom1)*fyi + a_Row(9,atom1)*fzi |
| 847 |
|
| 848 |
fxjj = a_Col(1,atom2)*fxj + a_Col(4,atom2)*fyj + a_Col(7,atom2)*fzj |
| 849 |
fyjj = a_Col(2,atom2)*fxj + a_Col(5,atom2)*fyj + a_Col(8,atom2)*fzj |
| 850 |
fzjj = a_Col(3,atom2)*fxj + a_Col(6,atom2)*fyj + a_Col(9,atom2)*fzj |
| 851 |
#else |
| 852 |
fxii = a(1,atom1)*fxi + a(4,atom1)*fyi + a(7,atom1)*fzi |
| 853 |
fyii = a(2,atom1)*fxi + a(5,atom1)*fyi + a(8,atom1)*fzi |
| 854 |
fzii = a(3,atom1)*fxi + a(6,atom1)*fyi + a(9,atom1)*fzi |
| 855 |
|
| 856 |
fxjj = a(1,atom2)*fxj + a(4,atom2)*fyj + a(7,atom2)*fzj |
| 857 |
fyjj = a(2,atom2)*fxj + a(5,atom2)*fyj + a(8,atom2)*fzj |
| 858 |
fzjj = a(3,atom2)*fxj + a(6,atom2)*fyj + a(9,atom2)*fzj |
| 859 |
#endif |
| 860 |
|
| 861 |
fxij = -fxii |
| 862 |
fyij = -fyii |
| 863 |
fzij = -fzii |
| 864 |
|
| 865 |
fxji = -fxjj |
| 866 |
fyji = -fyjj |
| 867 |
fzji = -fzjj |
| 868 |
|
| 869 |
fxradial = fxii + fxji |
| 870 |
fyradial = fyii + fyji |
| 871 |
fzradial = fzii + fzji |
| 872 |
|
| 873 |
#ifdef IS_MPI |
| 874 |
f_Row(1,atom1) = f_Row(1,atom1) + fxradial |
| 875 |
f_Row(2,atom1) = f_Row(2,atom1) + fyradial |
| 876 |
f_Row(3,atom1) = f_Row(3,atom1) + fzradial |
| 877 |
|
| 878 |
f_Col(1,atom2) = f_Col(1,atom2) - fxradial |
| 879 |
f_Col(2,atom2) = f_Col(2,atom2) - fyradial |
| 880 |
f_Col(3,atom2) = f_Col(3,atom2) - fzradial |
| 881 |
#else |
| 882 |
f(1,atom1) = f(1,atom1) + fxradial |
| 883 |
f(2,atom1) = f(2,atom1) + fyradial |
| 884 |
f(3,atom1) = f(3,atom1) + fzradial |
| 885 |
|
| 886 |
f(1,atom2) = f(1,atom2) - fxradial |
| 887 |
f(2,atom2) = f(2,atom2) - fyradial |
| 888 |
f(3,atom2) = f(3,atom2) - fzradial |
| 889 |
#endif |
| 890 |
|
| 891 |
#ifdef IS_MPI |
| 892 |
id1 = AtomRowToGlobal(atom1) |
| 893 |
id2 = AtomColToGlobal(atom2) |
| 894 |
#else |
| 895 |
id1 = atom1 |
| 896 |
id2 = atom2 |
| 897 |
#endif |
| 898 |
|
| 899 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
| 900 |
|
| 901 |
fpair(1) = fpair(1) + fxradial |
| 902 |
fpair(2) = fpair(2) + fyradial |
| 903 |
fpair(3) = fpair(3) + fzradial |
| 904 |
|
| 905 |
endif |
| 906 |
|
| 907 |
end subroutine do_shape_pair |
| 908 |
|
| 909 |
SUBROUTINE Associated_Legendre(x, l, m, lmax, plm, dlm) |
| 910 |
|
| 911 |
! Purpose: Compute the associated Legendre functions |
| 912 |
! Plm(x) and their derivatives Plm'(x) |
| 913 |
! Input : x --- Argument of Plm(x) |
| 914 |
! l --- Order of Plm(x), l = 0,1,2,...,n |
| 915 |
! m --- Degree of Plm(x), m = 0,1,2,...,N |
| 916 |
! lmax --- Physical dimension of PLM and DLM |
| 917 |
! Output: PLM(l,m) --- Plm(x) |
| 918 |
! DLM(l,m) --- Plm'(x) |
| 919 |
! |
| 920 |
! adapted from the routines in |
| 921 |
! COMPUTATION OF SPECIAL FUNCTIONS by Shanjie Zhang and Jianming Jin |
| 922 |
! ISBN 0-471-11963-6 |
| 923 |
! |
| 924 |
! The original Fortran77 codes can be found here: |
| 925 |
! http://iris-lee3.ece.uiuc.edu/~jjin/routines/routines.html |
| 926 |
|
| 927 |
real (kind=8), intent(in) :: x |
| 928 |
integer, intent(in) :: l, m, lmax |
| 929 |
real (kind=8), dimension(0:lmax,0:m), intent(out) :: PLM, DLM |
| 930 |
integer :: i, j, ls |
| 931 |
real (kind=8) :: xq, xs |
| 932 |
|
| 933 |
! zero out both arrays: |
| 934 |
DO I = 0, m |
| 935 |
DO J = 0, l |
| 936 |
PLM(J,I) = 0.0D0 |
| 937 |
DLM(J,I) = 0.0D0 |
| 938 |
end DO |
| 939 |
end DO |
| 940 |
|
| 941 |
! start with 0,0: |
| 942 |
PLM(0,0) = 1.0D0 |
| 943 |
|
| 944 |
! x = +/- 1 functions are easy: |
| 945 |
IF (abs(X).EQ.1.0D0) THEN |
| 946 |
DO I = 1, m |
| 947 |
PLM(0, I) = X**I |
| 948 |
DLM(0, I) = 0.5D0*I*(I+1.0D0)*X**(I+1) |
| 949 |
end DO |
| 950 |
DO J = 1, m |
| 951 |
DO I = 1, l |
| 952 |
IF (I.EQ.1) THEN |
| 953 |
DLM(I, J) = 1.0D+300 |
| 954 |
ELSE IF (I.EQ.2) THEN |
| 955 |
DLM(I, J) = -0.25D0*(J+2)*(J+1)*J*(J-1)*X**(J+1) |
| 956 |
ENDIF |
| 957 |
end DO |
| 958 |
end DO |
| 959 |
RETURN |
| 960 |
ENDIF |
| 961 |
|
| 962 |
LS = 1 |
| 963 |
IF (abs(X).GT.1.0D0) LS = -1 |
| 964 |
XQ = sqrt(LS*(1.0D0-X*X)) |
| 965 |
XS = LS*(1.0D0-X*X) |
| 966 |
|
| 967 |
DO I = 1, l |
| 968 |
PLM(I, I) = -LS*(2.0D0*I-1.0D0)*XQ*PLM(I-1, I-1) |
| 969 |
enddo |
| 970 |
|
| 971 |
DO I = 0, l |
| 972 |
PLM(I, I+1)=(2.0D0*I+1.0D0)*X*PLM(I, I) |
| 973 |
enddo |
| 974 |
|
| 975 |
DO I = 0, l |
| 976 |
DO J = I+2, m |
| 977 |
PLM(I, J)=((2.0D0*J-1.0D0)*X*PLM(I,J-1) - & |
| 978 |
(I+J-1.0D0)*PLM(I,J-2))/(J-I) |
| 979 |
end DO |
| 980 |
end DO |
| 981 |
|
| 982 |
DLM(0, 0)=0.0D0 |
| 983 |
|
| 984 |
DO J = 1, m |
| 985 |
DLM(0, J)=LS*J*(PLM(0,J-1)-X*PLM(0,J))/XS |
| 986 |
end DO |
| 987 |
|
| 988 |
DO I = 1, l |
| 989 |
DO J = I, m |
| 990 |
DLM(I,J) = LS*I*X*PLM(I, J)/XS + (J+I)*(J-I+1.0D0)/XQ*PLM(I-1, J) |
| 991 |
end DO |
| 992 |
end DO |
| 993 |
|
| 994 |
RETURN |
| 995 |
END SUBROUTINE Associated_Legendre |
| 996 |
|
| 997 |
|
| 998 |
subroutine Orthogonal_Polynomial(x, m, function_type, pl, dpl) |
| 999 |
|
| 1000 |
! Purpose: Compute orthogonal polynomials: Tn(x) or Un(x), |
| 1001 |
! or Ln(x) or Hn(x), and their derivatives |
| 1002 |
! Input : function_type --- Function code |
| 1003 |
! =1 for Chebyshev polynomial Tn(x) |
| 1004 |
! =2 for Chebyshev polynomial Un(x) |
| 1005 |
! =3 for Laguerre polynomial Ln(x) |
| 1006 |
! =4 for Hermite polynomial Hn(x) |
| 1007 |
! n --- Order of orthogonal polynomials |
| 1008 |
! x --- Argument of orthogonal polynomials |
| 1009 |
! Output: PL(n) --- Tn(x) or Un(x) or Ln(x) or Hn(x) |
| 1010 |
! DPL(n)--- Tn'(x) or Un'(x) or Ln'(x) or Hn'(x) |
| 1011 |
! |
| 1012 |
! adapted from the routines in |
| 1013 |
! COMPUTATION OF SPECIAL FUNCTIONS by Shanjie Zhang and Jianming Jin |
| 1014 |
! ISBN 0-471-11963-6 |
| 1015 |
! |
| 1016 |
! The original Fortran77 codes can be found here: |
| 1017 |
! http://iris-lee3.ece.uiuc.edu/~jjin/routines/routines.html |
| 1018 |
|
| 1019 |
real(kind=8), intent(in) :: x |
| 1020 |
integer, intent(in):: m |
| 1021 |
integer, intent(in):: function_type |
| 1022 |
real(kind=8), dimension(0:m), intent(inout) :: pl, dpl |
| 1023 |
|
| 1024 |
real(kind=8) :: a, b, c, y0, y1, dy0, dy1, yn, dyn |
| 1025 |
integer :: k |
| 1026 |
|
| 1027 |
A = 2.0D0 |
| 1028 |
B = 0.0D0 |
| 1029 |
C = 1.0D0 |
| 1030 |
Y0 = 1.0D0 |
| 1031 |
Y1 = 2.0D0*X |
| 1032 |
DY0 = 0.0D0 |
| 1033 |
DY1 = 2.0D0 |
| 1034 |
PL(0) = 1.0D0 |
| 1035 |
PL(1) = 2.0D0*X |
| 1036 |
DPL(0) = 0.0D0 |
| 1037 |
DPL(1) = 2.0D0 |
| 1038 |
IF (function_type.EQ.CHEBYSHEV_TN) THEN |
| 1039 |
Y1 = X |
| 1040 |
DY1 = 1.0D0 |
| 1041 |
PL(1) = X |
| 1042 |
DPL(1) = 1.0D0 |
| 1043 |
ELSE IF (function_type.EQ.LAGUERRE) THEN |
| 1044 |
Y1 = 1.0D0-X |
| 1045 |
DY1 = -1.0D0 |
| 1046 |
PL(1) = 1.0D0-X |
| 1047 |
DPL(1) = -1.0D0 |
| 1048 |
ENDIF |
| 1049 |
DO K = 2, m |
| 1050 |
IF (function_type.EQ.LAGUERRE) THEN |
| 1051 |
A = -1.0D0/K |
| 1052 |
B = 2.0D0+A |
| 1053 |
C = 1.0D0+A |
| 1054 |
ELSE IF (function_type.EQ.HERMITE) THEN |
| 1055 |
C = 2.0D0*(K-1.0D0) |
| 1056 |
ENDIF |
| 1057 |
YN = (A*X+B)*Y1-C*Y0 |
| 1058 |
DYN = A*Y1+(A*X+B)*DY1-C*DY0 |
| 1059 |
PL(K) = YN |
| 1060 |
DPL(K) = DYN |
| 1061 |
Y0 = Y1 |
| 1062 |
Y1 = YN |
| 1063 |
DY0 = DY1 |
| 1064 |
DY1 = DYN |
| 1065 |
end DO |
| 1066 |
RETURN |
| 1067 |
|
| 1068 |
end subroutine Orthogonal_Polynomial |
| 1069 |
|
| 1070 |
end module shapes |