| 1 |
gezelter |
1318 |
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 |
gezelter |
1321 |
|
| 28 |
|
|
type, private :: Shape |
| 29 |
|
|
integer :: atid |
| 30 |
|
|
integer :: nContactFuncs |
| 31 |
|
|
integer :: nRangeFuncs |
| 32 |
|
|
integer :: nStrengthFuncs |
| 33 |
|
|
integer :: bigL |
| 34 |
|
|
integer :: bigM |
| 35 |
|
|
integer, pointer, dimension(:) :: ContactFuncLValue => null() |
| 36 |
|
|
integer, pointer, dimension(:) :: ContactFuncMValue => null() |
| 37 |
|
|
integer, pointer, dimension(:) :: ContactFunctionType => null() |
| 38 |
|
|
real(kind=dp), pointer, dimension(:) :: ContactFuncCoefficient => null() |
| 39 |
|
|
integer, pointer, dimension(:) :: RangeFuncLValue => null() |
| 40 |
|
|
integer, pointer, dimension(:) :: RangeFuncMValue => null() |
| 41 |
|
|
integer, pointer, dimension(:) :: RangeFunctionType => null() |
| 42 |
|
|
real(kind=dp), pointer, dimension(:) :: RangeFuncCoefficient => null() |
| 43 |
|
|
integer, pointer, dimension(:) :: StrengthFuncLValue => null() |
| 44 |
|
|
integer, pointer, dimension(:) :: StrengthFuncMValue => null() |
| 45 |
|
|
integer, pointer, dimension(:) :: StrengthFunctionType => null() |
| 46 |
|
|
real(kind=dp), pointer, dimension(:) :: StrengthFuncCoefficient => null() |
| 47 |
|
|
logical :: isLJ |
| 48 |
|
|
real ( kind = dp ) :: epsilon |
| 49 |
|
|
real ( kind = dp ) :: sigma |
| 50 |
|
|
end type Shape |
| 51 |
|
|
|
| 52 |
|
|
type, private :: ShapeList |
| 53 |
|
|
integer :: n_shapes = 0 |
| 54 |
|
|
integer :: currentShape = 0 |
| 55 |
|
|
type (Shape), pointer :: Shapes(:) => null() |
| 56 |
|
|
integer, pointer :: atidToShape(:) => null() |
| 57 |
gezelter |
1318 |
end type ShapeList |
| 58 |
gezelter |
1321 |
|
| 59 |
|
|
type(ShapeList), save :: ShapeMap |
| 60 |
gezelter |
1318 |
|
| 61 |
|
|
integer :: lmax |
| 62 |
|
|
real (kind=dp), allocatable, dimension(:,:) :: plm_i, dlm_i, plm_j, dlm_j |
| 63 |
|
|
real (kind=dp), allocatable, dimension(:) :: tm_i, dtm_i, um_i, dum_i |
| 64 |
|
|
real (kind=dp), allocatable, dimension(:) :: tm_j, dtm_j, um_j, dum_j |
| 65 |
|
|
|
| 66 |
|
|
contains |
| 67 |
|
|
|
| 68 |
gezelter |
1321 |
subroutine newShapeType(nContactFuncs, ContactFuncLValue, & |
| 69 |
|
|
ContactFuncMValue, ContactFunctionType, ContactFuncCoefficient, & |
| 70 |
|
|
nRangeFuncs, RangeFuncLValue, RangeFuncMValue, RangeFunctionType, & |
| 71 |
|
|
RangeFuncCoefficient, nStrengthFuncs, StrengthFuncLValue, & |
| 72 |
|
|
StrengthFuncMValue, StrengthFunctionType, StrengthFuncCoefficient, & |
| 73 |
|
|
myAtid, status) |
| 74 |
|
|
|
| 75 |
|
|
integer :: nContactFuncs |
| 76 |
|
|
integer :: nRangeFuncs |
| 77 |
|
|
integer :: nStrengthFuncs |
| 78 |
|
|
integer :: shape_ident |
| 79 |
gezelter |
1318 |
integer :: status |
| 80 |
gezelter |
1321 |
integer :: myAtid |
| 81 |
|
|
integer :: bigL |
| 82 |
|
|
integer :: bigM |
| 83 |
|
|
integer :: j, me, nShapeTypes, nLJTypes, ntypes, current, alloc_stat |
| 84 |
|
|
integer, pointer :: MatchList(:) => null() |
| 85 |
|
|
|
| 86 |
|
|
integer, dimension(nContactFuncs) :: ContactFuncLValue |
| 87 |
|
|
integer, dimension(nContactFuncs) :: ContactFuncMValue |
| 88 |
|
|
integer, dimension(nContactFuncs) :: ContactFunctionType |
| 89 |
|
|
real(kind=dp), dimension(nContactFuncs) :: ContactFuncCoefficient |
| 90 |
|
|
integer, dimension(nRangeFuncs) :: RangeFuncLValue |
| 91 |
|
|
integer, dimension(nRangeFuncs) :: RangeFuncMValue |
| 92 |
|
|
integer, dimension(nRangeFuncs) :: RangeFunctionType |
| 93 |
|
|
real(kind=dp), dimension(nRangeFuncs) :: RangeFuncCoefficient |
| 94 |
|
|
integer, dimension(nStrengthFuncs) :: StrengthFuncLValue |
| 95 |
|
|
integer, dimension(nStrengthFuncs) :: StrengthFuncMValue |
| 96 |
|
|
integer, dimension(nStrengthFuncs) :: StrengthFunctionType |
| 97 |
|
|
real(kind=dp), dimension(nStrengthFuncs) :: StrengthFuncCoefficient |
| 98 |
|
|
|
| 99 |
|
|
status = 0 |
| 100 |
|
|
! check to see if this is the first time into this routine... |
| 101 |
|
|
if (.not.associated(ShapeMap%Shapes)) then |
| 102 |
|
|
|
| 103 |
|
|
call getMatchingElementList(atypes, "is_Shape", .true., & |
| 104 |
|
|
nShapeTypes, MatchList) |
| 105 |
|
|
|
| 106 |
|
|
call getMatchingElementList(atypes, "is_LJ", .true., & |
| 107 |
|
|
nLJTypes, MatchList) |
| 108 |
|
|
|
| 109 |
|
|
ShapeMap%n_shapes = nShapeTypes + nLJTypes |
| 110 |
|
|
|
| 111 |
|
|
allocate(ShapeMap%Shapes(nShapeTypes + nLJTypes)) |
| 112 |
|
|
|
| 113 |
|
|
ntypes = getSize(atypes) |
| 114 |
|
|
|
| 115 |
|
|
allocate(ShapeMap%atidToShape(ntypes)) |
| 116 |
|
|
end if |
| 117 |
|
|
|
| 118 |
|
|
ShapeMap%currentShape = ShapeMap%currentShape + 1 |
| 119 |
|
|
current = ShapeMap%currentShape |
| 120 |
|
|
|
| 121 |
|
|
call allocateShape(nContactFuncs, nRangeFuncs, nStrengthFuncs, & |
| 122 |
|
|
ShapeMap%Shapes(current), stat=alloc_stat) |
| 123 |
|
|
if (alloc_stat .ne. 0) then |
| 124 |
|
|
status = -1 |
| 125 |
|
|
return |
| 126 |
|
|
endif |
| 127 |
|
|
|
| 128 |
|
|
call getElementProperty(atypes, myAtid, "c_ident", me) |
| 129 |
|
|
ShapeMap%atidToShape(me) = current |
| 130 |
|
|
ShapeMap%Shapes(current)%atid = me |
| 131 |
|
|
ShapeMap%Shapes(current)%nContactFuncs = nContactFuncs |
| 132 |
|
|
ShapeMap%Shapes(current)%nRangeFuncs = nRangeFuncs |
| 133 |
|
|
ShapeMap%Shapes(current)%nStrengthFuncs = nStrengthFuncs |
| 134 |
|
|
ShapeMap%Shapes(current)%ContactFuncLValue = ContactFuncLValue |
| 135 |
|
|
ShapeMap%Shapes(current)%ContactFuncMValue = ContactFuncMValue |
| 136 |
|
|
ShapeMap%Shapes(current)%ContactFunctionType = ContactFunctionType |
| 137 |
|
|
ShapeMap%Shapes(current)%ContactFuncCoefficient = ContactFuncCoefficient |
| 138 |
|
|
ShapeMap%Shapes(current)%RangeFuncLValue = RangeFuncLValue |
| 139 |
|
|
ShapeMap%Shapes(current)%RangeFuncMValue = RangeFuncMValue |
| 140 |
|
|
ShapeMap%Shapes(current)%RangeFunctionType = RangeFunctionType |
| 141 |
|
|
ShapeMap%Shapes(current)%RangeFuncCoefficient = RangeFuncCoefficient |
| 142 |
|
|
ShapeMap%Shapes(current)%StrengthFuncLValue = StrengthFuncLValue |
| 143 |
|
|
ShapeMap%Shapes(current)%StrengthFuncMValue = StrengthFuncMValue |
| 144 |
|
|
ShapeMap%Shapes(current)%StrengthFunctionType = StrengthFunctionType |
| 145 |
|
|
ShapeMap%Shapes(current)%StrengthFuncCoefficient = StrengthFuncCoefficient |
| 146 |
|
|
|
| 147 |
|
|
bigL = -1 |
| 148 |
|
|
bigM = -1 |
| 149 |
|
|
|
| 150 |
|
|
do j = 1, ShapeMap%Shapes(current)%nContactFuncs |
| 151 |
|
|
if (ShapeMap%Shapes(current)%ContactFuncLValue(j) .gt. bigL) then |
| 152 |
|
|
bigL = ShapeMap%Shapes(current)%ContactFuncLValue(j) |
| 153 |
|
|
endif |
| 154 |
|
|
if (ShapeMap%Shapes(current)%ContactFuncMValue(j) .gt. bigM) then |
| 155 |
|
|
bigM = ShapeMap%Shapes(current)%ContactFuncMValue(j) |
| 156 |
|
|
endif |
| 157 |
|
|
enddo |
| 158 |
|
|
do j = 1, ShapeMap%Shapes(current)%nRangeFuncs |
| 159 |
|
|
if (ShapeMap%Shapes(current)%RangeFuncLValue(j) .gt. bigL) then |
| 160 |
|
|
bigL = ShapeMap%Shapes(current)%RangeFuncLValue(j) |
| 161 |
|
|
endif |
| 162 |
|
|
if (ShapeMap%Shapes(current)%RangeFuncMValue(j) .gt. bigM) then |
| 163 |
|
|
bigM = ShapeMap%Shapes(current)%RangeFuncMValue(j) |
| 164 |
|
|
endif |
| 165 |
|
|
enddo |
| 166 |
|
|
do j = 1, ShapeMap%Shapes(current)%nStrengthFuncs |
| 167 |
|
|
if (ShapeMap%Shapes(current)%StrengthFuncLValue(j) .gt. bigL) then |
| 168 |
|
|
bigL = ShapeMap%Shapes(current)%StrengthFuncLValue(j) |
| 169 |
|
|
endif |
| 170 |
|
|
if (ShapeMap%Shapes(current)%StrengthFuncMValue(j) .gt. bigM) then |
| 171 |
|
|
bigM = ShapeMap%Shapes(current)%StrengthFuncMValue(j) |
| 172 |
|
|
endif |
| 173 |
|
|
enddo |
| 174 |
|
|
|
| 175 |
|
|
ShapeMap%Shapes(current)%bigL = bigL |
| 176 |
|
|
ShapeMap%Shapes(current)%bigM = bigM |
| 177 |
|
|
|
| 178 |
|
|
end subroutine newShapeType |
| 179 |
|
|
|
| 180 |
|
|
subroutine allocateShape(nContactFuncs, nRangeFuncs, nStrengthFuncs, & |
| 181 |
|
|
myShape, stat) |
| 182 |
|
|
|
| 183 |
|
|
integer, intent(in) :: nContactFuncs, nRangeFuncs, nStrengthFuncs |
| 184 |
|
|
type(Shape), intent(inout) :: myShape |
| 185 |
|
|
integer, intent(out) :: stat |
| 186 |
|
|
integer :: alloc_stat |
| 187 |
|
|
|
| 188 |
|
|
if (associated(myShape%contactFuncLValue)) then |
| 189 |
|
|
deallocate(myShape%contactFuncLValue) |
| 190 |
|
|
endif |
| 191 |
|
|
allocate(myShape%contactFuncLValue(nContactFuncs), stat = alloc_stat) |
| 192 |
|
|
if (alloc_stat .ne. 0) then |
| 193 |
|
|
stat = -1 |
| 194 |
|
|
return |
| 195 |
|
|
endif |
| 196 |
|
|
if (associated(myShape%contactFuncMValue)) then |
| 197 |
|
|
deallocate(myShape%contactFuncMValue) |
| 198 |
|
|
endif |
| 199 |
|
|
allocate(myShape%contactFuncMValue(nContactFuncs), stat = alloc_stat) |
| 200 |
|
|
if (alloc_stat .ne. 0) then |
| 201 |
|
|
stat = -1 |
| 202 |
|
|
return |
| 203 |
|
|
endif |
| 204 |
|
|
if (associated(myShape%contactFunctionType)) then |
| 205 |
|
|
deallocate(myShape%contactFunctionType) |
| 206 |
|
|
endif |
| 207 |
|
|
allocate(myShape%contactFunctionType(nContactFuncs), stat = alloc_stat) |
| 208 |
|
|
if (alloc_stat .ne. 0) then |
| 209 |
|
|
stat = -1 |
| 210 |
|
|
return |
| 211 |
|
|
endif |
| 212 |
|
|
if (associated(myShape%contactFuncCoefficient)) then |
| 213 |
|
|
deallocate(myShape%contactFuncCoefficient) |
| 214 |
|
|
endif |
| 215 |
|
|
allocate(myShape%contactFuncCoefficient(nContactFuncs), stat = alloc_stat) |
| 216 |
|
|
if (alloc_stat .ne. 0) then |
| 217 |
|
|
stat = -1 |
| 218 |
|
|
return |
| 219 |
|
|
endif |
| 220 |
|
|
|
| 221 |
|
|
if (associated(myShape%rangeFuncLValue)) then |
| 222 |
|
|
deallocate(myShape%rangeFuncLValue) |
| 223 |
|
|
endif |
| 224 |
|
|
allocate(myShape%rangeFuncLValue(nRangeFuncs), stat = alloc_stat) |
| 225 |
|
|
if (alloc_stat .ne. 0) then |
| 226 |
|
|
stat = -1 |
| 227 |
|
|
return |
| 228 |
|
|
endif |
| 229 |
|
|
if (associated(myShape%rangeFuncMValue)) then |
| 230 |
|
|
deallocate(myShape%rangeFuncMValue) |
| 231 |
|
|
endif |
| 232 |
|
|
allocate(myShape%rangeFuncMValue(nRangeFuncs), stat = alloc_stat) |
| 233 |
|
|
if (alloc_stat .ne. 0) then |
| 234 |
|
|
stat = -1 |
| 235 |
|
|
return |
| 236 |
|
|
endif |
| 237 |
|
|
if (associated(myShape%rangeFunctionType)) then |
| 238 |
|
|
deallocate(myShape%rangeFunctionType) |
| 239 |
|
|
endif |
| 240 |
|
|
allocate(myShape%rangeFunctionType(nRangeFuncs), stat = alloc_stat) |
| 241 |
|
|
if (alloc_stat .ne. 0) then |
| 242 |
|
|
stat = -1 |
| 243 |
|
|
return |
| 244 |
|
|
endif |
| 245 |
|
|
if (associated(myShape%rangeFuncCoefficient)) then |
| 246 |
|
|
deallocate(myShape%rangeFuncCoefficient) |
| 247 |
|
|
endif |
| 248 |
|
|
allocate(myShape%rangeFuncCoefficient(nRangeFuncs), stat = alloc_stat) |
| 249 |
|
|
if (alloc_stat .ne. 0) then |
| 250 |
|
|
stat = -1 |
| 251 |
|
|
return |
| 252 |
|
|
endif |
| 253 |
|
|
|
| 254 |
|
|
if (associated(myShape%strengthFuncLValue)) then |
| 255 |
|
|
deallocate(myShape%strengthFuncLValue) |
| 256 |
|
|
endif |
| 257 |
|
|
allocate(myShape%strengthFuncLValue(nStrengthFuncs), stat = alloc_stat) |
| 258 |
|
|
if (alloc_stat .ne. 0) then |
| 259 |
|
|
stat = -1 |
| 260 |
|
|
return |
| 261 |
|
|
endif |
| 262 |
|
|
if (associated(myShape%strengthFuncMValue)) then |
| 263 |
|
|
deallocate(myShape%strengthFuncMValue) |
| 264 |
|
|
endif |
| 265 |
|
|
allocate(myShape%strengthFuncMValue(nStrengthFuncs), stat = alloc_stat) |
| 266 |
|
|
if (alloc_stat .ne. 0) then |
| 267 |
|
|
stat = -1 |
| 268 |
|
|
return |
| 269 |
|
|
endif |
| 270 |
|
|
if (associated(myShape%strengthFunctionType)) then |
| 271 |
|
|
deallocate(myShape%strengthFunctionType) |
| 272 |
|
|
endif |
| 273 |
|
|
allocate(myShape%strengthFunctionType(nStrengthFuncs), stat = alloc_stat) |
| 274 |
|
|
if (alloc_stat .ne. 0) then |
| 275 |
|
|
stat = -1 |
| 276 |
|
|
return |
| 277 |
|
|
endif |
| 278 |
|
|
if (associated(myShape%strengthFuncCoefficient)) then |
| 279 |
|
|
deallocate(myShape%strengthFuncCoefficient) |
| 280 |
|
|
endif |
| 281 |
|
|
allocate(myShape%strengthFuncCoefficient(nStrengthFuncs), stat=alloc_stat) |
| 282 |
|
|
if (alloc_stat .ne. 0) then |
| 283 |
|
|
stat = -1 |
| 284 |
|
|
return |
| 285 |
|
|
endif |
| 286 |
|
|
|
| 287 |
|
|
end subroutine allocateShape |
| 288 |
|
|
|
| 289 |
|
|
subroutine init_Shape_FF(status) |
| 290 |
|
|
integer :: status |
| 291 |
|
|
integer :: i, j, l, m, lm, function_type |
| 292 |
|
|
real(kind=dp) :: bigSigma, myBigSigma, thisSigma, coeff, Phunc, spi |
| 293 |
|
|
real(kind=dp) :: costheta, cpi, theta, Pi, phi, thisDP |
| 294 |
|
|
integer :: alloc_stat, iTheta, iPhi, nSteps, nAtypes, thisIP, current |
| 295 |
gezelter |
1318 |
logical :: thisProperty |
| 296 |
|
|
|
| 297 |
gezelter |
1321 |
Pi = 4.0d0 * datan(1.0d0) |
| 298 |
|
|
|
| 299 |
gezelter |
1318 |
status = 0 |
| 300 |
gezelter |
1321 |
if (ShapeMap%currentShape == 0) then |
| 301 |
|
|
call handleError("init_Shape_FF", "No members in ShapeMap") |
| 302 |
|
|
status = -1 |
| 303 |
|
|
return |
| 304 |
|
|
end if |
| 305 |
gezelter |
1318 |
|
| 306 |
gezelter |
1321 |
bigSigma = 0.0d0 |
| 307 |
|
|
do i = 1, ShapeMap%currentShape |
| 308 |
|
|
|
| 309 |
|
|
! Scan over theta and phi to |
| 310 |
|
|
! find the largest contact in any direction.... |
| 311 |
|
|
|
| 312 |
|
|
myBigSigma = 0.0d0 |
| 313 |
|
|
|
| 314 |
|
|
do iTheta = 0, nSteps |
| 315 |
|
|
theta = (Pi/2.0d0)*(dble(iTheta)/dble(nSteps)) |
| 316 |
|
|
costheta = cos(theta) |
| 317 |
|
|
|
| 318 |
|
|
call Associated_Legendre(costheta, ShapeMap%Shapes(i)%bigL, & |
| 319 |
|
|
ShapeMap%Shapes(i)%bigM, lmax, plm_i, dlm_i) |
| 320 |
|
|
|
| 321 |
|
|
do iPhi = 0, nSteps |
| 322 |
|
|
phi = -Pi + 2.0d0 * Pi * (dble(iPhi)/dble(nSteps)) |
| 323 |
|
|
cpi = cos(phi) |
| 324 |
|
|
spi = sin(phi) |
| 325 |
|
|
|
| 326 |
|
|
call Orthogonal_Polynomial(cpi, ShapeMap%Shapes(i)%bigM, & |
| 327 |
|
|
CHEBYSHEV_TN, tm_i, dtm_i) |
| 328 |
|
|
call Orthogonal_Polynomial(cpi, ShapeMap%Shapes(i)%bigM, & |
| 329 |
|
|
CHEBYSHEV_UN, um_i, dum_i) |
| 330 |
|
|
|
| 331 |
|
|
thisSigma = 0.0d0 |
| 332 |
|
|
|
| 333 |
|
|
do lm = 1, ShapeMap%Shapes(i)%nContactFuncs |
| 334 |
|
|
|
| 335 |
|
|
l = ShapeMap%Shapes(i)%ContactFuncLValue(lm) |
| 336 |
|
|
m = ShapeMap%Shapes(i)%ContactFuncMValue(lm) |
| 337 |
|
|
coeff = ShapeMap%Shapes(i)%ContactFuncCoefficient(lm) |
| 338 |
|
|
function_type = ShapeMap%Shapes(i)%ContactFunctionType(lm) |
| 339 |
|
|
|
| 340 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 341 |
|
|
Phunc = coeff * tm_i(m) |
| 342 |
|
|
else |
| 343 |
|
|
Phunc = coeff * spi * um_i(m-1) |
| 344 |
|
|
endif |
| 345 |
|
|
|
| 346 |
|
|
thisSigma = thisSigma + plm_i(l,m)*Phunc |
| 347 |
|
|
enddo |
| 348 |
|
|
|
| 349 |
|
|
if (thisSigma.gt.myBigSigma) myBigSigma = thisSigma |
| 350 |
|
|
enddo |
| 351 |
|
|
enddo |
| 352 |
|
|
|
| 353 |
|
|
if (myBigSigma.gt.bigSigma) bigSigma = myBigSigma |
| 354 |
|
|
enddo |
| 355 |
|
|
|
| 356 |
gezelter |
1318 |
nAtypes = getSize(atypes) |
| 357 |
|
|
|
| 358 |
|
|
if (nAtypes == 0) then |
| 359 |
|
|
status = -1 |
| 360 |
|
|
return |
| 361 |
|
|
end if |
| 362 |
|
|
|
| 363 |
|
|
do i = 1, nAtypes |
| 364 |
gezelter |
1321 |
|
| 365 |
|
|
call getElementProperty(atypes, i, "is_LJ", thisProperty) |
| 366 |
|
|
|
| 367 |
gezelter |
1318 |
if (thisProperty) then |
| 368 |
gezelter |
1321 |
|
| 369 |
|
|
ShapeMap%currentShape = ShapeMap%currentShape + 1 |
| 370 |
|
|
current = ShapeMap%currentShape |
| 371 |
gezelter |
1318 |
|
| 372 |
gezelter |
1321 |
call getElementProperty(atypes, i, "c_ident", thisIP) |
| 373 |
|
|
ShapeMap%atidToShape(thisIP) = current |
| 374 |
|
|
ShapeMap%Shapes(current)%atid = thisIP |
| 375 |
gezelter |
1318 |
|
| 376 |
gezelter |
1321 |
ShapeMap%Shapes(current)%isLJ = .true. |
| 377 |
gezelter |
1318 |
|
| 378 |
gezelter |
1321 |
call getElementProperty(atypes, i, "lj_epsilon", thisDP) |
| 379 |
|
|
ShapeMap%Shapes(current)%epsilon = thisDP |
| 380 |
gezelter |
1318 |
|
| 381 |
|
|
call getElementProperty(atypes, i, "lj_sigma", thisDP) |
| 382 |
gezelter |
1321 |
ShapeMap%Shapes(current)%sigma = thisDP |
| 383 |
|
|
if (thisDP .gt. bigSigma) bigSigma = thisDP |
| 384 |
|
|
|
| 385 |
gezelter |
1318 |
endif |
| 386 |
gezelter |
1321 |
|
| 387 |
gezelter |
1318 |
end do |
| 388 |
|
|
|
| 389 |
|
|
haveShapeMap = .true. |
| 390 |
gezelter |
1321 |
|
| 391 |
|
|
end subroutine init_Shape_FF |
| 392 |
|
|
|
| 393 |
gezelter |
1318 |
subroutine do_shape_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
| 394 |
|
|
pot, A, f, t, do_pot) |
| 395 |
|
|
|
| 396 |
|
|
integer, intent(in) :: atom1, atom2 |
| 397 |
|
|
real (kind=dp), intent(inout) :: rij, r2 |
| 398 |
|
|
real (kind=dp), dimension(3), intent(in) :: d |
| 399 |
|
|
real (kind=dp), dimension(3), intent(inout) :: fpair |
| 400 |
|
|
real (kind=dp) :: pot, vpair, sw |
| 401 |
|
|
real (kind=dp), dimension(9,nLocal) :: A |
| 402 |
|
|
real (kind=dp), dimension(3,nLocal) :: f |
| 403 |
|
|
real (kind=dp), dimension(3,nLocal) :: t |
| 404 |
|
|
logical, intent(in) :: do_pot |
| 405 |
|
|
|
| 406 |
|
|
real (kind=dp) :: r3, r5, rt2, rt3, rt5, rt6, rt11, rt12, rt126 |
| 407 |
gezelter |
1321 |
integer :: atid1, atid2, st1, st2 |
| 408 |
|
|
integer :: l, m, lm, id1, id2, localError, function_type |
| 409 |
gezelter |
1318 |
real (kind=dp) :: sigma_i, s_i, eps_i, sigma_j, s_j, eps_j |
| 410 |
|
|
real (kind=dp) :: coeff |
| 411 |
|
|
|
| 412 |
|
|
real (kind=dp) :: dsigmaidx, dsigmaidy, dsigmaidz |
| 413 |
|
|
real (kind=dp) :: dsigmaidux, dsigmaiduy, dsigmaiduz |
| 414 |
|
|
real (kind=dp) :: dsigmajdx, dsigmajdy, dsigmajdz |
| 415 |
|
|
real (kind=dp) :: dsigmajdux, dsigmajduy, dsigmajduz |
| 416 |
|
|
|
| 417 |
|
|
real (kind=dp) :: dsidx, dsidy, dsidz |
| 418 |
|
|
real (kind=dp) :: dsidux, dsiduy, dsiduz |
| 419 |
|
|
real (kind=dp) :: dsjdx, dsjdy, dsjdz |
| 420 |
|
|
real (kind=dp) :: dsjdux, dsjduy, dsjduz |
| 421 |
|
|
|
| 422 |
|
|
real (kind=dp) :: depsidx, depsidy, depsidz |
| 423 |
|
|
real (kind=dp) :: depsidux, depsiduy, depsiduz |
| 424 |
|
|
real (kind=dp) :: depsjdx, depsjdy, depsjdz |
| 425 |
|
|
real (kind=dp) :: depsjdux, depsjduy, depsjduz |
| 426 |
|
|
|
| 427 |
|
|
real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2 |
| 428 |
|
|
|
| 429 |
|
|
real (kind=dp) :: proji, proji3, projj, projj3 |
| 430 |
|
|
real (kind=dp) :: cti, ctj, cpi, cpj, spi, spj |
| 431 |
|
|
real (kind=dp) :: Phunc, sigma, s, eps, rtdenom, rt |
| 432 |
|
|
|
| 433 |
|
|
real (kind=dp) :: dctidx, dctidy, dctidz |
| 434 |
|
|
real (kind=dp) :: dctidux, dctiduy, dctiduz |
| 435 |
|
|
real (kind=dp) :: dctjdx, dctjdy, dctjdz |
| 436 |
|
|
real (kind=dp) :: dctjdux, dctjduy, dctjduz |
| 437 |
|
|
|
| 438 |
|
|
real (kind=dp) :: dcpidx, dcpidy, dcpidz |
| 439 |
|
|
real (kind=dp) :: dcpidux, dcpiduy, dcpiduz |
| 440 |
|
|
real (kind=dp) :: dcpjdx, dcpjdy, dcpjdz |
| 441 |
|
|
real (kind=dp) :: dcpjdux, dcpjduy, dcpjduz |
| 442 |
|
|
|
| 443 |
|
|
real (kind=dp) :: dspidx, dspidy, dspidz |
| 444 |
|
|
real (kind=dp) :: dspidux, dspiduy, dspiduz |
| 445 |
|
|
real (kind=dp) :: dspjdx, dspjdy, dspjdz |
| 446 |
|
|
real (kind=dp) :: dspjdux, dspjduy, dspjduz |
| 447 |
|
|
|
| 448 |
|
|
real (kind=dp) :: dPhuncdX, dPhuncdY, dPhuncdZ |
| 449 |
|
|
real (kind=dp) :: dPhuncdUx, dPhuncdUy, dPhuncdUz |
| 450 |
|
|
|
| 451 |
|
|
real (kind=dp) :: dsigmadxi, dsigmadyi, dsigmadzi |
| 452 |
|
|
real (kind=dp) :: dsigmaduxi, dsigmaduyi, dsigmaduzi |
| 453 |
|
|
real (kind=dp) :: dsigmadxj, dsigmadyj, dsigmadzj |
| 454 |
|
|
real (kind=dp) :: dsigmaduxj, dsigmaduyj, dsigmaduzj |
| 455 |
|
|
|
| 456 |
|
|
real (kind=dp) :: dsdxi, dsdyi, dsdzi |
| 457 |
|
|
real (kind=dp) :: dsduxi, dsduyi, dsduzi |
| 458 |
|
|
real (kind=dp) :: dsdxj, dsdyj, dsdzj |
| 459 |
|
|
real (kind=dp) :: dsduxj, dsduyj, dsduzj |
| 460 |
|
|
|
| 461 |
|
|
real (kind=dp) :: depsdxi, depsdyi, depsdzi |
| 462 |
|
|
real (kind=dp) :: depsduxi, depsduyi, depsduzi |
| 463 |
|
|
real (kind=dp) :: depsdxj, depsdyj, depsdzj |
| 464 |
|
|
real (kind=dp) :: depsduxj, depsduyj, depsduzj |
| 465 |
|
|
|
| 466 |
|
|
real (kind=dp) :: drtdxi, drtdyi, drtdzi |
| 467 |
|
|
real (kind=dp) :: drtduxi, drtduyi, drtduzi |
| 468 |
|
|
real (kind=dp) :: drtdxj, drtdyj, drtdzj |
| 469 |
|
|
real (kind=dp) :: drtduxj, drtduyj, drtduzj |
| 470 |
|
|
|
| 471 |
|
|
real (kind=dp) :: drdxi, drdyi, drdzi |
| 472 |
|
|
real (kind=dp) :: drduxi, drduyi, drduzi |
| 473 |
|
|
real (kind=dp) :: drdxj, drdyj, drdzj |
| 474 |
|
|
real (kind=dp) :: drduxj, drduyj, drduzj |
| 475 |
|
|
|
| 476 |
|
|
real (kind=dp) :: dvdxi, dvdyi, dvdzi |
| 477 |
|
|
real (kind=dp) :: dvduxi, dvduyi, dvduzi |
| 478 |
|
|
real (kind=dp) :: dvdxj, dvdyj, dvdzj |
| 479 |
|
|
real (kind=dp) :: dvduxj, dvduyj, dvduzj |
| 480 |
|
|
|
| 481 |
|
|
real (kind=dp) :: fxi, fyi, fzi, fxj, fyj, fzj |
| 482 |
|
|
real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj |
| 483 |
|
|
real (kind=dp) :: fxii, fyii, fzii, fxij, fyij, fzij |
| 484 |
|
|
real (kind=dp) :: fxji, fyji, fzji, fxjj, fyjj, fzjj |
| 485 |
|
|
real (kind=dp) :: fxradial, fyradial, fzradial |
| 486 |
|
|
|
| 487 |
|
|
if (.not.haveShapeMap) then |
| 488 |
gezelter |
1321 |
call handleError("calc_shape", "NO SHAPEMAP!!!!") |
| 489 |
|
|
return |
| 490 |
gezelter |
1318 |
endif |
| 491 |
|
|
|
| 492 |
|
|
!! We assume that the rotation matrices have already been calculated |
| 493 |
|
|
!! and placed in the A array. |
| 494 |
|
|
|
| 495 |
|
|
r3 = r2*rij |
| 496 |
|
|
r5 = r3*r2 |
| 497 |
|
|
|
| 498 |
|
|
drdxi = -d(1) / rij |
| 499 |
|
|
drdyi = -d(2) / rij |
| 500 |
|
|
drdzi = -d(3) / rij |
| 501 |
|
|
|
| 502 |
|
|
drdxj = d(1) / rij |
| 503 |
|
|
drdyj = d(2) / rij |
| 504 |
|
|
drdzj = d(3) / rij |
| 505 |
|
|
|
| 506 |
gezelter |
1321 |
! find the atom type id (atid) for each atom: |
| 507 |
gezelter |
1318 |
#ifdef IS_MPI |
| 508 |
gezelter |
1321 |
atid1 = atid_Row(atom1) |
| 509 |
|
|
atid2 = atid_Col(atom2) |
| 510 |
gezelter |
1318 |
#else |
| 511 |
gezelter |
1321 |
atid1 = atid(atom1) |
| 512 |
|
|
atid2 = atid(atom2) |
| 513 |
gezelter |
1318 |
#endif |
| 514 |
|
|
|
| 515 |
gezelter |
1321 |
! use the atid to find the shape type (st) for each atom: |
| 516 |
|
|
|
| 517 |
|
|
st1 = ShapeMap%atidToShape(atid1) |
| 518 |
|
|
st2 = ShapeMap%atidToShape(atid2) |
| 519 |
|
|
|
| 520 |
|
|
if (ShapeMap%Shapes(st1)%isLJ) then |
| 521 |
|
|
sigma_i = ShapeMap%Shapes(st1)%sigma |
| 522 |
|
|
s_i = ShapeMap%Shapes(st1)%sigma |
| 523 |
|
|
eps_i = ShapeMap%Shapes(st1)%epsilon |
| 524 |
gezelter |
1318 |
dsigmaidx = 0.0d0 |
| 525 |
|
|
dsigmaidy = 0.0d0 |
| 526 |
|
|
dsigmaidz = 0.0d0 |
| 527 |
|
|
dsigmaidux = 0.0d0 |
| 528 |
|
|
dsigmaiduy = 0.0d0 |
| 529 |
|
|
dsigmaiduz = 0.0d0 |
| 530 |
|
|
dsidx = 0.0d0 |
| 531 |
|
|
dsidy = 0.0d0 |
| 532 |
|
|
dsidz = 0.0d0 |
| 533 |
|
|
dsidux = 0.0d0 |
| 534 |
|
|
dsiduy = 0.0d0 |
| 535 |
|
|
dsiduz = 0.0d0 |
| 536 |
|
|
depsidx = 0.0d0 |
| 537 |
|
|
depsidy = 0.0d0 |
| 538 |
|
|
depsidz = 0.0d0 |
| 539 |
|
|
depsidux = 0.0d0 |
| 540 |
|
|
depsiduy = 0.0d0 |
| 541 |
|
|
depsiduz = 0.0d0 |
| 542 |
|
|
else |
| 543 |
|
|
|
| 544 |
|
|
#ifdef IS_MPI |
| 545 |
|
|
! rotate the inter-particle separation into the two different |
| 546 |
|
|
! body-fixed coordinate systems: |
| 547 |
|
|
|
| 548 |
|
|
xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3) |
| 549 |
|
|
yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3) |
| 550 |
|
|
zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3) |
| 551 |
|
|
|
| 552 |
|
|
#else |
| 553 |
|
|
! rotate the inter-particle separation into the two different |
| 554 |
|
|
! body-fixed coordinate systems: |
| 555 |
|
|
|
| 556 |
|
|
xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3) |
| 557 |
|
|
yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3) |
| 558 |
|
|
zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3) |
| 559 |
|
|
|
| 560 |
|
|
#endif |
| 561 |
|
|
|
| 562 |
|
|
xi2 = xi*xi |
| 563 |
|
|
yi2 = yi*yi |
| 564 |
|
|
zi2 = zi*zi |
| 565 |
|
|
|
| 566 |
|
|
proji = sqrt(xi2 + yi2) |
| 567 |
|
|
proji3 = proji*proji*proji |
| 568 |
|
|
|
| 569 |
|
|
cti = zi / rij |
| 570 |
|
|
dctidx = - zi * xi / r3 |
| 571 |
|
|
dctidy = - zi * yi / r3 |
| 572 |
|
|
dctidz = 1.0d0 / rij - zi2 / r3 |
| 573 |
|
|
dctidux = yi / rij |
| 574 |
|
|
dctiduy = -xi / rij |
| 575 |
|
|
dctiduz = 0.0d0 |
| 576 |
|
|
|
| 577 |
|
|
cpi = xi / proji |
| 578 |
|
|
dcpidx = 1.0d0 / proji - xi2 / proji3 |
| 579 |
|
|
dcpidy = - xi * yi / proji3 |
| 580 |
|
|
dcpidz = 0.0d0 |
| 581 |
|
|
dcpidux = xi * yi * zi / proji3 |
| 582 |
|
|
dcpiduy = -zi * (1.0d0 / proji - xi2 / proji3) |
| 583 |
|
|
dcpiduz = -yi * (1.0d0 / proji - xi2 / proji3) - (xi2 * yi / proji3) |
| 584 |
|
|
|
| 585 |
|
|
spi = yi / proji |
| 586 |
|
|
dspidx = - xi * yi / proji3 |
| 587 |
|
|
dspidy = 1.0d0 / proji - yi2 / proji3 |
| 588 |
|
|
dspidz = 0.0d0 |
| 589 |
|
|
dspidux = -zi * (1.0d0 / proji - yi2 / proji3) |
| 590 |
|
|
dspiduy = xi * yi * zi / proji3 |
| 591 |
|
|
dspiduz = xi * (1.0d0 / proji - yi2 / proji3) + (xi * yi2 / proji3) |
| 592 |
|
|
|
| 593 |
gezelter |
1321 |
call Associated_Legendre(cti, ShapeMap%Shapes(st1)%bigL, & |
| 594 |
|
|
ShapeMap%Shapes(st1)%bigM, lmax, plm_i, dlm_i) |
| 595 |
gezelter |
1318 |
|
| 596 |
gezelter |
1321 |
call Orthogonal_Polynomial(cpi, ShapeMap%Shapes(st1)%bigM, & |
| 597 |
|
|
CHEBYSHEV_TN, tm_i, dtm_i) |
| 598 |
|
|
call Orthogonal_Polynomial(cpi, ShapeMap%Shapes(st1)%bigM, & |
| 599 |
|
|
CHEBYSHEV_UN, um_i, dum_i) |
| 600 |
gezelter |
1318 |
|
| 601 |
|
|
sigma_i = 0.0d0 |
| 602 |
|
|
s_i = 0.0d0 |
| 603 |
|
|
eps_i = 0.0d0 |
| 604 |
|
|
dsigmaidx = 0.0d0 |
| 605 |
|
|
dsigmaidy = 0.0d0 |
| 606 |
|
|
dsigmaidz = 0.0d0 |
| 607 |
|
|
dsigmaidux = 0.0d0 |
| 608 |
|
|
dsigmaiduy = 0.0d0 |
| 609 |
|
|
dsigmaiduz = 0.0d0 |
| 610 |
|
|
dsidx = 0.0d0 |
| 611 |
|
|
dsidy = 0.0d0 |
| 612 |
|
|
dsidz = 0.0d0 |
| 613 |
|
|
dsidux = 0.0d0 |
| 614 |
|
|
dsiduy = 0.0d0 |
| 615 |
|
|
dsiduz = 0.0d0 |
| 616 |
|
|
depsidx = 0.0d0 |
| 617 |
|
|
depsidy = 0.0d0 |
| 618 |
|
|
depsidz = 0.0d0 |
| 619 |
|
|
depsidux = 0.0d0 |
| 620 |
|
|
depsiduy = 0.0d0 |
| 621 |
|
|
depsiduz = 0.0d0 |
| 622 |
|
|
|
| 623 |
gezelter |
1321 |
do lm = 1, ShapeMap%Shapes(st1)%nContactFuncs |
| 624 |
|
|
l = ShapeMap%Shapes(st1)%ContactFuncLValue(lm) |
| 625 |
|
|
m = ShapeMap%Shapes(st1)%ContactFuncMValue(lm) |
| 626 |
|
|
coeff = ShapeMap%Shapes(st1)%ContactFuncCoefficient(lm) |
| 627 |
|
|
function_type = ShapeMap%Shapes(st1)%ContactFunctionType(lm) |
| 628 |
gezelter |
1318 |
|
| 629 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 630 |
|
|
Phunc = coeff * tm_i(m) |
| 631 |
|
|
dPhuncdX = coeff * dtm_i(m) * dcpidx |
| 632 |
|
|
dPhuncdY = coeff * dtm_i(m) * dcpidy |
| 633 |
|
|
dPhuncdZ = coeff * dtm_i(m) * dcpidz |
| 634 |
|
|
dPhuncdUz = coeff * dtm_i(m) * dcpidux |
| 635 |
|
|
dPhuncdUy = coeff * dtm_i(m) * dcpiduy |
| 636 |
|
|
dPhuncdUz = coeff * dtm_i(m) * dcpiduz |
| 637 |
|
|
else |
| 638 |
|
|
Phunc = coeff * spi * um_i(m-1) |
| 639 |
|
|
dPhuncdX = coeff * (spi * dum_i(m-1) * dcpidx + dspidx *um_i(m-1)) |
| 640 |
|
|
dPhuncdY = coeff * (spi * dum_i(m-1) * dcpidy + dspidy *um_i(m-1)) |
| 641 |
|
|
dPhuncdZ = coeff * (spi * dum_i(m-1) * dcpidz + dspidz *um_i(m-1)) |
| 642 |
|
|
dPhuncdUx = coeff*(spi * dum_i(m-1)*dcpidux + dspidux *um_i(m-1)) |
| 643 |
|
|
dPhuncdUy = coeff*(spi * dum_i(m-1)*dcpiduy + dspiduy *um_i(m-1)) |
| 644 |
|
|
dPhuncdUz = coeff*(spi * dum_i(m-1)*dcpiduz + dspiduz *um_i(m-1)) |
| 645 |
|
|
endif |
| 646 |
|
|
|
| 647 |
|
|
sigma_i = sigma_i + plm_i(l,m)*Phunc |
| 648 |
|
|
|
| 649 |
|
|
dsigmaidx = dsigmaidx + plm_i(l,m)*dPhuncdX + & |
| 650 |
|
|
Phunc * dlm_i(l,m) * dctidx |
| 651 |
|
|
dsigmaidy = dsigmaidy + plm_i(l,m)*dPhuncdY + & |
| 652 |
|
|
Phunc * dlm_i(l,m) * dctidy |
| 653 |
|
|
dsigmaidz = dsigmaidz + plm_i(l,m)*dPhuncdZ + & |
| 654 |
|
|
Phunc * dlm_i(l,m) * dctidz |
| 655 |
|
|
|
| 656 |
|
|
dsigmaidux = dsigmaidux + plm_i(l,m)* dPhuncdUx + & |
| 657 |
|
|
Phunc * dlm_i(l,m) * dctidux |
| 658 |
|
|
dsigmaiduy = dsigmaiduy + plm_i(l,m)* dPhuncdUy + & |
| 659 |
|
|
Phunc * dlm_i(l,m) * dctiduy |
| 660 |
|
|
dsigmaiduz = dsigmaiduz + plm_i(l,m)* dPhuncdUz + & |
| 661 |
|
|
Phunc * dlm_i(l,m) * dctiduz |
| 662 |
|
|
|
| 663 |
|
|
end do |
| 664 |
|
|
|
| 665 |
gezelter |
1321 |
do lm = 1, ShapeMap%Shapes(st1)%nRangeFuncs |
| 666 |
|
|
l = ShapeMap%Shapes(st1)%RangeFuncLValue(lm) |
| 667 |
|
|
m = ShapeMap%Shapes(st1)%RangeFuncMValue(lm) |
| 668 |
|
|
coeff = ShapeMap%Shapes(st1)%RangeFuncCoefficient(lm) |
| 669 |
|
|
function_type = ShapeMap%Shapes(st1)%RangeFunctionType(lm) |
| 670 |
gezelter |
1318 |
|
| 671 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 672 |
|
|
Phunc = coeff * tm_i(m) |
| 673 |
|
|
dPhuncdX = coeff * dtm_i(m) * dcpidx |
| 674 |
|
|
dPhuncdY = coeff * dtm_i(m) * dcpidy |
| 675 |
|
|
dPhuncdZ = coeff * dtm_i(m) * dcpidz |
| 676 |
|
|
dPhuncdUz = coeff * dtm_i(m) * dcpidux |
| 677 |
|
|
dPhuncdUy = coeff * dtm_i(m) * dcpiduy |
| 678 |
|
|
dPhuncdUz = coeff * dtm_i(m) * dcpiduz |
| 679 |
|
|
else |
| 680 |
|
|
Phunc = coeff * spi * um_i(m-1) |
| 681 |
|
|
dPhuncdX = coeff * (spi * dum_i(m-1) * dcpidx + dspidx *um_i(m-1)) |
| 682 |
|
|
dPhuncdY = coeff * (spi * dum_i(m-1) * dcpidy + dspidy *um_i(m-1)) |
| 683 |
|
|
dPhuncdZ = coeff * (spi * dum_i(m-1) * dcpidz + dspidz *um_i(m-1)) |
| 684 |
|
|
dPhuncdUx = coeff*(spi * dum_i(m-1)*dcpidux + dspidux *um_i(m-1)) |
| 685 |
|
|
dPhuncdUy = coeff*(spi * dum_i(m-1)*dcpiduy + dspiduy *um_i(m-1)) |
| 686 |
|
|
dPhuncdUz = coeff*(spi * dum_i(m-1)*dcpiduz + dspiduz *um_i(m-1)) |
| 687 |
|
|
endif |
| 688 |
|
|
|
| 689 |
|
|
s_i = s_i + plm_i(l,m)*Phunc |
| 690 |
|
|
|
| 691 |
|
|
dsidx = dsidx + plm_i(l,m)*dPhuncdX + & |
| 692 |
|
|
Phunc * dlm_i(l,m) * dctidx |
| 693 |
|
|
dsidy = dsidy + plm_i(l,m)*dPhuncdY + & |
| 694 |
|
|
Phunc * dlm_i(l,m) * dctidy |
| 695 |
|
|
dsidz = dsidz + plm_i(l,m)*dPhuncdZ + & |
| 696 |
|
|
Phunc * dlm_i(l,m) * dctidz |
| 697 |
|
|
|
| 698 |
|
|
dsidux = dsidux + plm_i(l,m)* dPhuncdUx + & |
| 699 |
|
|
Phunc * dlm_i(l,m) * dctidux |
| 700 |
|
|
dsiduy = dsiduy + plm_i(l,m)* dPhuncdUy + & |
| 701 |
|
|
Phunc * dlm_i(l,m) * dctiduy |
| 702 |
|
|
dsiduz = dsiduz + plm_i(l,m)* dPhuncdUz + & |
| 703 |
|
|
Phunc * dlm_i(l,m) * dctiduz |
| 704 |
|
|
|
| 705 |
|
|
end do |
| 706 |
|
|
|
| 707 |
gezelter |
1321 |
do lm = 1, ShapeMap%Shapes(st1)%nStrengthFuncs |
| 708 |
|
|
l = ShapeMap%Shapes(st1)%StrengthFuncLValue(lm) |
| 709 |
|
|
m = ShapeMap%Shapes(st1)%StrengthFuncMValue(lm) |
| 710 |
|
|
coeff = ShapeMap%Shapes(st1)%StrengthFuncCoefficient(lm) |
| 711 |
|
|
function_type = ShapeMap%Shapes(st1)%StrengthFunctionType(lm) |
| 712 |
gezelter |
1318 |
|
| 713 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 714 |
|
|
Phunc = coeff * tm_i(m) |
| 715 |
|
|
dPhuncdX = coeff * dtm_i(m) * dcpidx |
| 716 |
|
|
dPhuncdY = coeff * dtm_i(m) * dcpidy |
| 717 |
|
|
dPhuncdZ = coeff * dtm_i(m) * dcpidz |
| 718 |
|
|
dPhuncdUz = coeff * dtm_i(m) * dcpidux |
| 719 |
|
|
dPhuncdUy = coeff * dtm_i(m) * dcpiduy |
| 720 |
|
|
dPhuncdUz = coeff * dtm_i(m) * dcpiduz |
| 721 |
|
|
else |
| 722 |
|
|
Phunc = coeff * spi * um_i(m-1) |
| 723 |
|
|
dPhuncdX = coeff * (spi * dum_i(m-1) * dcpidx + dspidx *um_i(m-1)) |
| 724 |
|
|
dPhuncdY = coeff * (spi * dum_i(m-1) * dcpidy + dspidy *um_i(m-1)) |
| 725 |
|
|
dPhuncdZ = coeff * (spi * dum_i(m-1) * dcpidz + dspidz *um_i(m-1)) |
| 726 |
|
|
dPhuncdUx = coeff*(spi * dum_i(m-1)*dcpidux + dspidux *um_i(m-1)) |
| 727 |
|
|
dPhuncdUy = coeff*(spi * dum_i(m-1)*dcpiduy + dspiduy *um_i(m-1)) |
| 728 |
|
|
dPhuncdUz = coeff*(spi * dum_i(m-1)*dcpiduz + dspiduz *um_i(m-1)) |
| 729 |
|
|
endif |
| 730 |
|
|
|
| 731 |
|
|
eps_i = eps_i + plm_i(l,m)*Phunc |
| 732 |
|
|
|
| 733 |
|
|
depsidx = depsidx + plm_i(l,m)*dPhuncdX + & |
| 734 |
|
|
Phunc * dlm_i(l,m) * dctidx |
| 735 |
|
|
depsidy = depsidy + plm_i(l,m)*dPhuncdY + & |
| 736 |
|
|
Phunc * dlm_i(l,m) * dctidy |
| 737 |
|
|
depsidz = depsidz + plm_i(l,m)*dPhuncdZ + & |
| 738 |
|
|
Phunc * dlm_i(l,m) * dctidz |
| 739 |
|
|
|
| 740 |
|
|
depsidux = depsidux + plm_i(l,m)* dPhuncdUx + & |
| 741 |
|
|
Phunc * dlm_i(l,m) * dctidux |
| 742 |
|
|
depsiduy = depsiduy + plm_i(l,m)* dPhuncdUy + & |
| 743 |
|
|
Phunc * dlm_i(l,m) * dctiduy |
| 744 |
|
|
depsiduz = depsiduz + plm_i(l,m)* dPhuncdUz + & |
| 745 |
|
|
Phunc * dlm_i(l,m) * dctiduz |
| 746 |
|
|
|
| 747 |
|
|
end do |
| 748 |
|
|
|
| 749 |
|
|
endif |
| 750 |
|
|
|
| 751 |
|
|
! now do j: |
| 752 |
|
|
|
| 753 |
gezelter |
1321 |
if (ShapeMap%Shapes(st2)%isLJ) then |
| 754 |
|
|
sigma_j = ShapeMap%Shapes(st2)%sigma |
| 755 |
|
|
s_j = ShapeMap%Shapes(st2)%sigma |
| 756 |
|
|
eps_j = ShapeMap%Shapes(st2)%epsilon |
| 757 |
gezelter |
1318 |
dsigmajdx = 0.0d0 |
| 758 |
|
|
dsigmajdy = 0.0d0 |
| 759 |
|
|
dsigmajdz = 0.0d0 |
| 760 |
|
|
dsigmajdux = 0.0d0 |
| 761 |
|
|
dsigmajduy = 0.0d0 |
| 762 |
|
|
dsigmajduz = 0.0d0 |
| 763 |
|
|
dsjdx = 0.0d0 |
| 764 |
|
|
dsjdy = 0.0d0 |
| 765 |
|
|
dsjdz = 0.0d0 |
| 766 |
|
|
dsjdux = 0.0d0 |
| 767 |
|
|
dsjduy = 0.0d0 |
| 768 |
|
|
dsjduz = 0.0d0 |
| 769 |
|
|
depsjdx = 0.0d0 |
| 770 |
|
|
depsjdy = 0.0d0 |
| 771 |
|
|
depsjdz = 0.0d0 |
| 772 |
|
|
depsjdux = 0.0d0 |
| 773 |
|
|
depsjduy = 0.0d0 |
| 774 |
|
|
depsjduz = 0.0d0 |
| 775 |
|
|
else |
| 776 |
|
|
|
| 777 |
|
|
#ifdef IS_MPI |
| 778 |
|
|
! rotate the inter-particle separation into the two different |
| 779 |
|
|
! body-fixed coordinate systems: |
| 780 |
|
|
! negative sign because this is the vector from j to i: |
| 781 |
|
|
|
| 782 |
|
|
xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3)) |
| 783 |
|
|
yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3)) |
| 784 |
|
|
zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3)) |
| 785 |
|
|
#else |
| 786 |
|
|
! rotate the inter-particle separation into the two different |
| 787 |
|
|
! body-fixed coordinate systems: |
| 788 |
|
|
! negative sign because this is the vector from j to i: |
| 789 |
|
|
|
| 790 |
|
|
xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3)) |
| 791 |
|
|
yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3)) |
| 792 |
|
|
zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3)) |
| 793 |
|
|
#endif |
| 794 |
|
|
|
| 795 |
|
|
xj2 = xj*xj |
| 796 |
|
|
yj2 = yj*yj |
| 797 |
|
|
zj2 = zj*zj |
| 798 |
|
|
|
| 799 |
|
|
projj = sqrt(xj2 + yj2) |
| 800 |
|
|
projj3 = projj*projj*projj |
| 801 |
|
|
|
| 802 |
|
|
ctj = zj / rij |
| 803 |
|
|
dctjdx = - zj * xj / r3 |
| 804 |
|
|
dctjdy = - zj * yj / r3 |
| 805 |
|
|
dctjdz = 1.0d0 / rij - zj2 / r3 |
| 806 |
|
|
dctjdux = yj / rij |
| 807 |
|
|
dctjduy = -xj / rij |
| 808 |
|
|
dctjduz = 0.0d0 |
| 809 |
|
|
|
| 810 |
|
|
cpj = xj / projj |
| 811 |
|
|
dcpjdx = 1.0d0 / projj - xj2 / projj3 |
| 812 |
|
|
dcpjdy = - xj * yj / projj3 |
| 813 |
|
|
dcpjdz = 0.0d0 |
| 814 |
|
|
dcpjdux = xj * yj * zj / projj3 |
| 815 |
|
|
dcpjduy = -zj * (1.0d0 / projj - xj2 / projj3) |
| 816 |
|
|
dcpjduz = -yj * (1.0d0 / projj - xj2 / projj3) - (xj2 * yj / projj3) |
| 817 |
|
|
|
| 818 |
|
|
spj = yj / projj |
| 819 |
|
|
dspjdx = - xj * yj / projj3 |
| 820 |
|
|
dspjdy = 1.0d0 / projj - yj2 / projj3 |
| 821 |
|
|
dspjdz = 0.0d0 |
| 822 |
|
|
dspjdux = -zj * (1.0d0 / projj - yj2 / projj3) |
| 823 |
|
|
dspjduy = xj * yj * zj / projj3 |
| 824 |
|
|
dspjduz = xj * (1.0d0 / projj - yi2 / projj3) + (xj * yj2 / projj3) |
| 825 |
|
|
|
| 826 |
gezelter |
1321 |
call Associated_Legendre(ctj, ShapeMap%Shapes(st2)%bigL, & |
| 827 |
|
|
ShapeMap%Shapes(st2)%bigM, lmax, plm_j, dlm_j) |
| 828 |
gezelter |
1318 |
|
| 829 |
gezelter |
1321 |
call Orthogonal_Polynomial(cpj, ShapeMap%Shapes(st2)%bigM, & |
| 830 |
|
|
CHEBYSHEV_TN, tm_j, dtm_j) |
| 831 |
|
|
call Orthogonal_Polynomial(cpj, ShapeMap%Shapes(st2)%bigM, & |
| 832 |
|
|
CHEBYSHEV_UN, um_j, dum_j) |
| 833 |
gezelter |
1318 |
|
| 834 |
|
|
sigma_j = 0.0d0 |
| 835 |
|
|
s_j = 0.0d0 |
| 836 |
|
|
eps_j = 0.0d0 |
| 837 |
|
|
dsigmajdx = 0.0d0 |
| 838 |
|
|
dsigmajdy = 0.0d0 |
| 839 |
|
|
dsigmajdz = 0.0d0 |
| 840 |
|
|
dsigmajdux = 0.0d0 |
| 841 |
|
|
dsigmajduy = 0.0d0 |
| 842 |
|
|
dsigmajduz = 0.0d0 |
| 843 |
|
|
dsjdx = 0.0d0 |
| 844 |
|
|
dsjdy = 0.0d0 |
| 845 |
|
|
dsjdz = 0.0d0 |
| 846 |
|
|
dsjdux = 0.0d0 |
| 847 |
|
|
dsjduy = 0.0d0 |
| 848 |
|
|
dsjduz = 0.0d0 |
| 849 |
|
|
depsjdx = 0.0d0 |
| 850 |
|
|
depsjdy = 0.0d0 |
| 851 |
|
|
depsjdz = 0.0d0 |
| 852 |
|
|
depsjdux = 0.0d0 |
| 853 |
|
|
depsjduy = 0.0d0 |
| 854 |
|
|
depsjduz = 0.0d0 |
| 855 |
|
|
|
| 856 |
gezelter |
1321 |
do lm = 1, ShapeMap%Shapes(st2)%nContactFuncs |
| 857 |
|
|
l = ShapeMap%Shapes(st2)%ContactFuncLValue(lm) |
| 858 |
|
|
m = ShapeMap%Shapes(st2)%ContactFuncMValue(lm) |
| 859 |
|
|
coeff = ShapeMap%Shapes(st2)%ContactFuncCoefficient(lm) |
| 860 |
|
|
function_type = ShapeMap%Shapes(st2)%ContactFunctionType(lm) |
| 861 |
gezelter |
1318 |
|
| 862 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 863 |
|
|
Phunc = coeff * tm_j(m) |
| 864 |
|
|
dPhuncdX = coeff * dtm_j(m) * dcpjdx |
| 865 |
|
|
dPhuncdY = coeff * dtm_j(m) * dcpjdy |
| 866 |
|
|
dPhuncdZ = coeff * dtm_j(m) * dcpjdz |
| 867 |
|
|
dPhuncdUz = coeff * dtm_j(m) * dcpjdux |
| 868 |
|
|
dPhuncdUy = coeff * dtm_j(m) * dcpjduy |
| 869 |
|
|
dPhuncdUz = coeff * dtm_j(m) * dcpjduz |
| 870 |
|
|
else |
| 871 |
|
|
Phunc = coeff * spj * um_j(m-1) |
| 872 |
|
|
dPhuncdX = coeff * (spj * dum_j(m-1) * dcpjdx + dspjdx *um_j(m-1)) |
| 873 |
|
|
dPhuncdY = coeff * (spj * dum_j(m-1) * dcpjdy + dspjdy *um_j(m-1)) |
| 874 |
|
|
dPhuncdZ = coeff * (spj * dum_j(m-1) * dcpjdz + dspjdz *um_j(m-1)) |
| 875 |
|
|
dPhuncdUx = coeff*(spj * dum_j(m-1)*dcpjdux + dspjdux *um_j(m-1)) |
| 876 |
|
|
dPhuncdUy = coeff*(spj * dum_j(m-1)*dcpjduy + dspjduy *um_j(m-1)) |
| 877 |
|
|
dPhuncdUz = coeff*(spj * dum_j(m-1)*dcpjduz + dspjduz *um_j(m-1)) |
| 878 |
|
|
endif |
| 879 |
|
|
|
| 880 |
|
|
sigma_j = sigma_j + plm_j(l,m)*Phunc |
| 881 |
|
|
|
| 882 |
|
|
dsigmajdx = dsigmajdx + plm_j(l,m)*dPhuncdX + & |
| 883 |
|
|
Phunc * dlm_j(l,m) * dctjdx |
| 884 |
|
|
dsigmajdy = dsigmajdy + plm_j(l,m)*dPhuncdY + & |
| 885 |
|
|
Phunc * dlm_j(l,m) * dctjdy |
| 886 |
|
|
dsigmajdz = dsigmajdz + plm_j(l,m)*dPhuncdZ + & |
| 887 |
|
|
Phunc * dlm_j(l,m) * dctjdz |
| 888 |
|
|
|
| 889 |
|
|
dsigmajdux = dsigmajdux + plm_j(l,m)* dPhuncdUx + & |
| 890 |
|
|
Phunc * dlm_j(l,m) * dctjdux |
| 891 |
|
|
dsigmajduy = dsigmajduy + plm_j(l,m)* dPhuncdUy + & |
| 892 |
|
|
Phunc * dlm_j(l,m) * dctjduy |
| 893 |
|
|
dsigmajduz = dsigmajduz + plm_j(l,m)* dPhuncdUz + & |
| 894 |
|
|
Phunc * dlm_j(l,m) * dctjduz |
| 895 |
|
|
|
| 896 |
|
|
end do |
| 897 |
|
|
|
| 898 |
gezelter |
1321 |
do lm = 1, ShapeMap%Shapes(st2)%nRangeFuncs |
| 899 |
|
|
l = ShapeMap%Shapes(st2)%RangeFuncLValue(lm) |
| 900 |
|
|
m = ShapeMap%Shapes(st2)%RangeFuncMValue(lm) |
| 901 |
|
|
coeff = ShapeMap%Shapes(st2)%RangeFuncCoefficient(lm) |
| 902 |
|
|
function_type = ShapeMap%Shapes(st2)%RangeFunctionType(lm) |
| 903 |
gezelter |
1318 |
|
| 904 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 905 |
|
|
Phunc = coeff * tm_j(m) |
| 906 |
|
|
dPhuncdX = coeff * dtm_j(m) * dcpjdx |
| 907 |
|
|
dPhuncdY = coeff * dtm_j(m) * dcpjdy |
| 908 |
|
|
dPhuncdZ = coeff * dtm_j(m) * dcpjdz |
| 909 |
|
|
dPhuncdUz = coeff * dtm_j(m) * dcpjdux |
| 910 |
|
|
dPhuncdUy = coeff * dtm_j(m) * dcpjduy |
| 911 |
|
|
dPhuncdUz = coeff * dtm_j(m) * dcpjduz |
| 912 |
|
|
else |
| 913 |
|
|
Phunc = coeff * spj * um_j(m-1) |
| 914 |
|
|
dPhuncdX = coeff * (spj * dum_j(m-1) * dcpjdx + dspjdx *um_j(m-1)) |
| 915 |
|
|
dPhuncdY = coeff * (spj * dum_j(m-1) * dcpjdy + dspjdy *um_j(m-1)) |
| 916 |
|
|
dPhuncdZ = coeff * (spj * dum_j(m-1) * dcpjdz + dspjdz *um_j(m-1)) |
| 917 |
|
|
dPhuncdUx = coeff*(spj * dum_j(m-1)*dcpjdux + dspjdux *um_j(m-1)) |
| 918 |
|
|
dPhuncdUy = coeff*(spj * dum_j(m-1)*dcpjduy + dspjduy *um_j(m-1)) |
| 919 |
|
|
dPhuncdUz = coeff*(spj * dum_j(m-1)*dcpjduz + dspjduz *um_j(m-1)) |
| 920 |
|
|
endif |
| 921 |
|
|
|
| 922 |
|
|
s_j = s_j + plm_j(l,m)*Phunc |
| 923 |
|
|
|
| 924 |
|
|
dsjdx = dsjdx + plm_j(l,m)*dPhuncdX + & |
| 925 |
|
|
Phunc * dlm_j(l,m) * dctjdx |
| 926 |
|
|
dsjdy = dsjdy + plm_j(l,m)*dPhuncdY + & |
| 927 |
|
|
Phunc * dlm_j(l,m) * dctjdy |
| 928 |
|
|
dsjdz = dsjdz + plm_j(l,m)*dPhuncdZ + & |
| 929 |
|
|
Phunc * dlm_j(l,m) * dctjdz |
| 930 |
|
|
|
| 931 |
|
|
dsjdux = dsjdux + plm_j(l,m)* dPhuncdUx + & |
| 932 |
|
|
Phunc * dlm_j(l,m) * dctjdux |
| 933 |
|
|
dsjduy = dsjduy + plm_j(l,m)* dPhuncdUy + & |
| 934 |
|
|
Phunc * dlm_j(l,m) * dctjduy |
| 935 |
|
|
dsjduz = dsjduz + plm_j(l,m)* dPhuncdUz + & |
| 936 |
|
|
Phunc * dlm_j(l,m) * dctjduz |
| 937 |
|
|
|
| 938 |
|
|
end do |
| 939 |
|
|
|
| 940 |
gezelter |
1321 |
do lm = 1, ShapeMap%Shapes(st2)%nStrengthFuncs |
| 941 |
|
|
l = ShapeMap%Shapes(st2)%StrengthFuncLValue(lm) |
| 942 |
|
|
m = ShapeMap%Shapes(st2)%StrengthFuncMValue(lm) |
| 943 |
|
|
coeff = ShapeMap%Shapes(st2)%StrengthFuncCoefficient(lm) |
| 944 |
|
|
function_type = ShapeMap%Shapes(st2)%StrengthFunctionType(lm) |
| 945 |
gezelter |
1318 |
|
| 946 |
|
|
if ((function_type .eq. SH_COS).or.(m.eq.0)) then |
| 947 |
|
|
Phunc = coeff * tm_j(m) |
| 948 |
|
|
dPhuncdX = coeff * dtm_j(m) * dcpjdx |
| 949 |
|
|
dPhuncdY = coeff * dtm_j(m) * dcpjdy |
| 950 |
|
|
dPhuncdZ = coeff * dtm_j(m) * dcpjdz |
| 951 |
|
|
dPhuncdUz = coeff * dtm_j(m) * dcpjdux |
| 952 |
|
|
dPhuncdUy = coeff * dtm_j(m) * dcpjduy |
| 953 |
|
|
dPhuncdUz = coeff * dtm_j(m) * dcpjduz |
| 954 |
|
|
else |
| 955 |
|
|
Phunc = coeff * spj * um_j(m-1) |
| 956 |
|
|
dPhuncdX = coeff * (spj * dum_j(m-1) * dcpjdx + dspjdx *um_j(m-1)) |
| 957 |
|
|
dPhuncdY = coeff * (spj * dum_j(m-1) * dcpjdy + dspjdy *um_j(m-1)) |
| 958 |
|
|
dPhuncdZ = coeff * (spj * dum_j(m-1) * dcpjdz + dspjdz *um_j(m-1)) |
| 959 |
|
|
dPhuncdUx = coeff*(spj * dum_j(m-1)*dcpjdux + dspjdux *um_j(m-1)) |
| 960 |
|
|
dPhuncdUy = coeff*(spj * dum_j(m-1)*dcpjduy + dspjduy *um_j(m-1)) |
| 961 |
|
|
dPhuncdUz = coeff*(spj * dum_j(m-1)*dcpjduz + dspjduz *um_j(m-1)) |
| 962 |
|
|
endif |
| 963 |
|
|
|
| 964 |
|
|
eps_j = eps_j + plm_j(l,m)*Phunc |
| 965 |
|
|
|
| 966 |
|
|
depsjdx = depsjdx + plm_j(l,m)*dPhuncdX + & |
| 967 |
|
|
Phunc * dlm_j(l,m) * dctjdx |
| 968 |
|
|
depsjdy = depsjdy + plm_j(l,m)*dPhuncdY + & |
| 969 |
|
|
Phunc * dlm_j(l,m) * dctjdy |
| 970 |
|
|
depsjdz = depsjdz + plm_j(l,m)*dPhuncdZ + & |
| 971 |
|
|
Phunc * dlm_j(l,m) * dctjdz |
| 972 |
|
|
|
| 973 |
|
|
depsjdux = depsjdux + plm_j(l,m)* dPhuncdUx + & |
| 974 |
|
|
Phunc * dlm_j(l,m) * dctjdux |
| 975 |
|
|
depsjduy = depsjduy + plm_j(l,m)* dPhuncdUy + & |
| 976 |
|
|
Phunc * dlm_j(l,m) * dctjduy |
| 977 |
|
|
depsjduz = depsjduz + plm_j(l,m)* dPhuncdUz + & |
| 978 |
|
|
Phunc * dlm_j(l,m) * dctjduz |
| 979 |
|
|
|
| 980 |
|
|
end do |
| 981 |
|
|
|
| 982 |
|
|
endif |
| 983 |
|
|
|
| 984 |
|
|
! phew, now let's assemble the potential energy: |
| 985 |
|
|
|
| 986 |
|
|
sigma = 0.5*(sigma_i + sigma_j) |
| 987 |
|
|
|
| 988 |
|
|
dsigmadxi = 0.5*dsigmaidx |
| 989 |
|
|
dsigmadyi = 0.5*dsigmaidy |
| 990 |
|
|
dsigmadzi = 0.5*dsigmaidz |
| 991 |
|
|
dsigmaduxi = 0.5*dsigmaidux |
| 992 |
|
|
dsigmaduyi = 0.5*dsigmaiduy |
| 993 |
|
|
dsigmaduzi = 0.5*dsigmaiduz |
| 994 |
|
|
|
| 995 |
|
|
dsigmadxj = 0.5*dsigmajdx |
| 996 |
|
|
dsigmadyj = 0.5*dsigmajdy |
| 997 |
|
|
dsigmadzj = 0.5*dsigmajdz |
| 998 |
|
|
dsigmaduxj = 0.5*dsigmajdux |
| 999 |
|
|
dsigmaduyj = 0.5*dsigmajduy |
| 1000 |
|
|
dsigmaduzj = 0.5*dsigmajduz |
| 1001 |
|
|
|
| 1002 |
|
|
s = 0.5*(s_i + s_j) |
| 1003 |
|
|
|
| 1004 |
|
|
dsdxi = 0.5*dsidx |
| 1005 |
|
|
dsdyi = 0.5*dsidy |
| 1006 |
|
|
dsdzi = 0.5*dsidz |
| 1007 |
|
|
dsduxi = 0.5*dsidux |
| 1008 |
|
|
dsduyi = 0.5*dsiduy |
| 1009 |
|
|
dsduzi = 0.5*dsiduz |
| 1010 |
|
|
|
| 1011 |
|
|
dsdxj = 0.5*dsjdx |
| 1012 |
|
|
dsdyj = 0.5*dsjdy |
| 1013 |
|
|
dsdzj = 0.5*dsjdz |
| 1014 |
|
|
dsduxj = 0.5*dsjdux |
| 1015 |
|
|
dsduyj = 0.5*dsjduy |
| 1016 |
|
|
dsduzj = 0.5*dsjduz |
| 1017 |
|
|
|
| 1018 |
|
|
eps = sqrt(eps_i * eps_j) |
| 1019 |
|
|
|
| 1020 |
|
|
depsdxi = eps_j * depsidx / (2.0d0 * eps) |
| 1021 |
|
|
depsdyi = eps_j * depsidy / (2.0d0 * eps) |
| 1022 |
|
|
depsdzi = eps_j * depsidz / (2.0d0 * eps) |
| 1023 |
|
|
depsduxi = eps_j * depsidux / (2.0d0 * eps) |
| 1024 |
|
|
depsduyi = eps_j * depsiduy / (2.0d0 * eps) |
| 1025 |
|
|
depsduzi = eps_j * depsiduz / (2.0d0 * eps) |
| 1026 |
|
|
|
| 1027 |
|
|
depsdxj = eps_i * depsjdx / (2.0d0 * eps) |
| 1028 |
|
|
depsdyj = eps_i * depsjdy / (2.0d0 * eps) |
| 1029 |
|
|
depsdzj = eps_i * depsjdz / (2.0d0 * eps) |
| 1030 |
|
|
depsduxj = eps_i * depsjdux / (2.0d0 * eps) |
| 1031 |
|
|
depsduyj = eps_i * depsjduy / (2.0d0 * eps) |
| 1032 |
|
|
depsduzj = eps_i * depsjduz / (2.0d0 * eps) |
| 1033 |
|
|
|
| 1034 |
|
|
rtdenom = rij-sigma+s |
| 1035 |
|
|
rt = s / rtdenom |
| 1036 |
|
|
|
| 1037 |
|
|
drtdxi = (dsdxi + rt * (drdxi - dsigmadxi + dsdxi)) / rtdenom |
| 1038 |
|
|
drtdyi = (dsdyi + rt * (drdyi - dsigmadyi + dsdyi)) / rtdenom |
| 1039 |
|
|
drtdzi = (dsdzi + rt * (drdzi - dsigmadzi + dsdzi)) / rtdenom |
| 1040 |
|
|
drtduxi = (dsduxi + rt * (drduxi - dsigmaduxi + dsduxi)) / rtdenom |
| 1041 |
|
|
drtduyi = (dsduyi + rt * (drduyi - dsigmaduyi + dsduyi)) / rtdenom |
| 1042 |
|
|
drtduzi = (dsduzi + rt * (drduzi - dsigmaduzi + dsduzi)) / rtdenom |
| 1043 |
|
|
drtdxj = (dsdxj + rt * (drdxj - dsigmadxj + dsdxj)) / rtdenom |
| 1044 |
|
|
drtdyj = (dsdyj + rt * (drdyj - dsigmadyj + dsdyj)) / rtdenom |
| 1045 |
|
|
drtdzj = (dsdzj + rt * (drdzj - dsigmadzj + dsdzj)) / rtdenom |
| 1046 |
|
|
drtduxj = (dsduxj + rt * (drduxj - dsigmaduxj + dsduxj)) / rtdenom |
| 1047 |
|
|
drtduyj = (dsduyj + rt * (drduyj - dsigmaduyj + dsduyj)) / rtdenom |
| 1048 |
|
|
drtduzj = (dsduzj + rt * (drduzj - dsigmaduzj + dsduzj)) / rtdenom |
| 1049 |
|
|
|
| 1050 |
|
|
rt2 = rt*rt |
| 1051 |
|
|
rt3 = rt2*rt |
| 1052 |
|
|
rt5 = rt2*rt3 |
| 1053 |
|
|
rt6 = rt3*rt3 |
| 1054 |
|
|
rt11 = rt5*rt6 |
| 1055 |
|
|
rt12 = rt6*rt6 |
| 1056 |
|
|
rt126 = rt12 - rt6 |
| 1057 |
|
|
|
| 1058 |
|
|
if (do_pot) then |
| 1059 |
|
|
#ifdef IS_MPI |
| 1060 |
|
|
pot_row(atom1) = pot_row(atom1) + 2.0d0*eps*rt126*sw |
| 1061 |
|
|
pot_col(atom2) = pot_col(atom2) + 2.0d0*eps*rt126*sw |
| 1062 |
|
|
#else |
| 1063 |
|
|
pot = pot + 4.0d0*eps*rt126*sw |
| 1064 |
|
|
#endif |
| 1065 |
|
|
endif |
| 1066 |
|
|
|
| 1067 |
|
|
dvdxi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdxi + 4.0d0*depsdxi*rt126 |
| 1068 |
|
|
dvdyi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdyi + 4.0d0*depsdyi*rt126 |
| 1069 |
|
|
dvdzi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdzi + 4.0d0*depsdzi*rt126 |
| 1070 |
|
|
dvduxi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduxi + 4.0d0*depsduxi*rt126 |
| 1071 |
|
|
dvduyi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduyi + 4.0d0*depsduyi*rt126 |
| 1072 |
|
|
dvduzi = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduzi + 4.0d0*depsduzi*rt126 |
| 1073 |
|
|
|
| 1074 |
|
|
dvdxj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdxj + 4.0d0*depsdxj*rt126 |
| 1075 |
|
|
dvdyj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdyj + 4.0d0*depsdyj*rt126 |
| 1076 |
|
|
dvdzj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtdzj + 4.0d0*depsdzj*rt126 |
| 1077 |
|
|
dvduxj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduxj + 4.0d0*depsduxj*rt126 |
| 1078 |
|
|
dvduyj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduyj + 4.0d0*depsduyj*rt126 |
| 1079 |
|
|
dvduzj = 24.0d0*eps*(2.0d0*rt11 - rt5)*drtduzj + 4.0d0*depsduzj*rt126 |
| 1080 |
|
|
|
| 1081 |
|
|
! do the torques first since they are easy: |
| 1082 |
|
|
! remember that these are still in the body fixed axes |
| 1083 |
|
|
|
| 1084 |
|
|
txi = dvduxi * sw |
| 1085 |
|
|
tyi = dvduyi * sw |
| 1086 |
|
|
tzi = dvduzi * sw |
| 1087 |
|
|
|
| 1088 |
|
|
txj = dvduxj * sw |
| 1089 |
|
|
tyj = dvduyj * sw |
| 1090 |
|
|
tzj = dvduzj * sw |
| 1091 |
|
|
|
| 1092 |
|
|
! go back to lab frame using transpose of rotation matrix: |
| 1093 |
|
|
|
| 1094 |
|
|
#ifdef IS_MPI |
| 1095 |
|
|
t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + & |
| 1096 |
|
|
a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi |
| 1097 |
|
|
t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + & |
| 1098 |
|
|
a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi |
| 1099 |
|
|
t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + & |
| 1100 |
|
|
a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi |
| 1101 |
|
|
|
| 1102 |
|
|
t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + & |
| 1103 |
|
|
a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj |
| 1104 |
|
|
t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + & |
| 1105 |
|
|
a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj |
| 1106 |
|
|
t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + & |
| 1107 |
|
|
a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj |
| 1108 |
|
|
#else |
| 1109 |
|
|
t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi |
| 1110 |
|
|
t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi |
| 1111 |
|
|
t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi |
| 1112 |
|
|
|
| 1113 |
|
|
t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj |
| 1114 |
|
|
t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj |
| 1115 |
|
|
t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj |
| 1116 |
|
|
#endif |
| 1117 |
|
|
! Now, on to the forces: |
| 1118 |
|
|
|
| 1119 |
|
|
! first rotate the i terms back into the lab frame: |
| 1120 |
|
|
|
| 1121 |
|
|
fxi = dvdxi * sw |
| 1122 |
|
|
fyi = dvdyi * sw |
| 1123 |
|
|
fzi = dvdzi * sw |
| 1124 |
|
|
|
| 1125 |
|
|
fxj = dvdxj * sw |
| 1126 |
|
|
fyj = dvdyj * sw |
| 1127 |
|
|
fzj = dvdzj * sw |
| 1128 |
|
|
|
| 1129 |
|
|
#ifdef IS_MPI |
| 1130 |
|
|
fxii = a_Row(1,atom1)*fxi + a_Row(4,atom1)*fyi + a_Row(7,atom1)*fzi |
| 1131 |
|
|
fyii = a_Row(2,atom1)*fxi + a_Row(5,atom1)*fyi + a_Row(8,atom1)*fzi |
| 1132 |
|
|
fzii = a_Row(3,atom1)*fxi + a_Row(6,atom1)*fyi + a_Row(9,atom1)*fzi |
| 1133 |
|
|
|
| 1134 |
|
|
fxjj = a_Col(1,atom2)*fxj + a_Col(4,atom2)*fyj + a_Col(7,atom2)*fzj |
| 1135 |
|
|
fyjj = a_Col(2,atom2)*fxj + a_Col(5,atom2)*fyj + a_Col(8,atom2)*fzj |
| 1136 |
|
|
fzjj = a_Col(3,atom2)*fxj + a_Col(6,atom2)*fyj + a_Col(9,atom2)*fzj |
| 1137 |
|
|
#else |
| 1138 |
|
|
fxii = a(1,atom1)*fxi + a(4,atom1)*fyi + a(7,atom1)*fzi |
| 1139 |
|
|
fyii = a(2,atom1)*fxi + a(5,atom1)*fyi + a(8,atom1)*fzi |
| 1140 |
|
|
fzii = a(3,atom1)*fxi + a(6,atom1)*fyi + a(9,atom1)*fzi |
| 1141 |
|
|
|
| 1142 |
|
|
fxjj = a(1,atom2)*fxj + a(4,atom2)*fyj + a(7,atom2)*fzj |
| 1143 |
|
|
fyjj = a(2,atom2)*fxj + a(5,atom2)*fyj + a(8,atom2)*fzj |
| 1144 |
|
|
fzjj = a(3,atom2)*fxj + a(6,atom2)*fyj + a(9,atom2)*fzj |
| 1145 |
|
|
#endif |
| 1146 |
|
|
|
| 1147 |
|
|
fxij = -fxii |
| 1148 |
|
|
fyij = -fyii |
| 1149 |
|
|
fzij = -fzii |
| 1150 |
|
|
|
| 1151 |
|
|
fxji = -fxjj |
| 1152 |
|
|
fyji = -fyjj |
| 1153 |
|
|
fzji = -fzjj |
| 1154 |
|
|
|
| 1155 |
|
|
fxradial = fxii + fxji |
| 1156 |
|
|
fyradial = fyii + fyji |
| 1157 |
|
|
fzradial = fzii + fzji |
| 1158 |
|
|
|
| 1159 |
|
|
#ifdef IS_MPI |
| 1160 |
|
|
f_Row(1,atom1) = f_Row(1,atom1) + fxradial |
| 1161 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + fyradial |
| 1162 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + fzradial |
| 1163 |
|
|
|
| 1164 |
|
|
f_Col(1,atom2) = f_Col(1,atom2) - fxradial |
| 1165 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - fyradial |
| 1166 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - fzradial |
| 1167 |
|
|
#else |
| 1168 |
|
|
f(1,atom1) = f(1,atom1) + fxradial |
| 1169 |
|
|
f(2,atom1) = f(2,atom1) + fyradial |
| 1170 |
|
|
f(3,atom1) = f(3,atom1) + fzradial |
| 1171 |
|
|
|
| 1172 |
|
|
f(1,atom2) = f(1,atom2) - fxradial |
| 1173 |
|
|
f(2,atom2) = f(2,atom2) - fyradial |
| 1174 |
|
|
f(3,atom2) = f(3,atom2) - fzradial |
| 1175 |
|
|
#endif |
| 1176 |
|
|
|
| 1177 |
|
|
#ifdef IS_MPI |
| 1178 |
|
|
id1 = AtomRowToGlobal(atom1) |
| 1179 |
|
|
id2 = AtomColToGlobal(atom2) |
| 1180 |
|
|
#else |
| 1181 |
|
|
id1 = atom1 |
| 1182 |
|
|
id2 = atom2 |
| 1183 |
|
|
#endif |
| 1184 |
|
|
|
| 1185 |
|
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
| 1186 |
|
|
|
| 1187 |
|
|
fpair(1) = fpair(1) + fxradial |
| 1188 |
|
|
fpair(2) = fpair(2) + fyradial |
| 1189 |
|
|
fpair(3) = fpair(3) + fzradial |
| 1190 |
|
|
|
| 1191 |
|
|
endif |
| 1192 |
|
|
|
| 1193 |
|
|
end subroutine do_shape_pair |
| 1194 |
|
|
|
| 1195 |
|
|
SUBROUTINE Associated_Legendre(x, l, m, lmax, plm, dlm) |
| 1196 |
|
|
|
| 1197 |
|
|
! Purpose: Compute the associated Legendre functions |
| 1198 |
|
|
! Plm(x) and their derivatives Plm'(x) |
| 1199 |
|
|
! Input : x --- Argument of Plm(x) |
| 1200 |
|
|
! l --- Order of Plm(x), l = 0,1,2,...,n |
| 1201 |
|
|
! m --- Degree of Plm(x), m = 0,1,2,...,N |
| 1202 |
|
|
! lmax --- Physical dimension of PLM and DLM |
| 1203 |
|
|
! Output: PLM(l,m) --- Plm(x) |
| 1204 |
|
|
! DLM(l,m) --- Plm'(x) |
| 1205 |
|
|
! |
| 1206 |
|
|
! adapted from the routines in |
| 1207 |
|
|
! COMPUTATION OF SPECIAL FUNCTIONS by Shanjie Zhang and Jianming Jin |
| 1208 |
|
|
! ISBN 0-471-11963-6 |
| 1209 |
|
|
! |
| 1210 |
|
|
! The original Fortran77 codes can be found here: |
| 1211 |
|
|
! http://iris-lee3.ece.uiuc.edu/~jjin/routines/routines.html |
| 1212 |
|
|
|
| 1213 |
|
|
real (kind=8), intent(in) :: x |
| 1214 |
|
|
integer, intent(in) :: l, m, lmax |
| 1215 |
|
|
real (kind=8), dimension(0:lmax,0:m), intent(out) :: PLM, DLM |
| 1216 |
|
|
integer :: i, j, ls |
| 1217 |
|
|
real (kind=8) :: xq, xs |
| 1218 |
|
|
|
| 1219 |
|
|
! zero out both arrays: |
| 1220 |
|
|
DO I = 0, m |
| 1221 |
|
|
DO J = 0, l |
| 1222 |
|
|
PLM(J,I) = 0.0D0 |
| 1223 |
|
|
DLM(J,I) = 0.0D0 |
| 1224 |
|
|
end DO |
| 1225 |
|
|
end DO |
| 1226 |
|
|
|
| 1227 |
|
|
! start with 0,0: |
| 1228 |
|
|
PLM(0,0) = 1.0D0 |
| 1229 |
|
|
|
| 1230 |
|
|
! x = +/- 1 functions are easy: |
| 1231 |
|
|
IF (abs(X).EQ.1.0D0) THEN |
| 1232 |
|
|
DO I = 1, m |
| 1233 |
|
|
PLM(0, I) = X**I |
| 1234 |
|
|
DLM(0, I) = 0.5D0*I*(I+1.0D0)*X**(I+1) |
| 1235 |
|
|
end DO |
| 1236 |
|
|
DO J = 1, m |
| 1237 |
|
|
DO I = 1, l |
| 1238 |
|
|
IF (I.EQ.1) THEN |
| 1239 |
|
|
DLM(I, J) = 1.0D+300 |
| 1240 |
|
|
ELSE IF (I.EQ.2) THEN |
| 1241 |
|
|
DLM(I, J) = -0.25D0*(J+2)*(J+1)*J*(J-1)*X**(J+1) |
| 1242 |
|
|
ENDIF |
| 1243 |
|
|
end DO |
| 1244 |
|
|
end DO |
| 1245 |
|
|
RETURN |
| 1246 |
|
|
ENDIF |
| 1247 |
|
|
|
| 1248 |
|
|
LS = 1 |
| 1249 |
|
|
IF (abs(X).GT.1.0D0) LS = -1 |
| 1250 |
|
|
XQ = sqrt(LS*(1.0D0-X*X)) |
| 1251 |
|
|
XS = LS*(1.0D0-X*X) |
| 1252 |
|
|
|
| 1253 |
|
|
DO I = 1, l |
| 1254 |
|
|
PLM(I, I) = -LS*(2.0D0*I-1.0D0)*XQ*PLM(I-1, I-1) |
| 1255 |
|
|
enddo |
| 1256 |
|
|
|
| 1257 |
|
|
DO I = 0, l |
| 1258 |
|
|
PLM(I, I+1)=(2.0D0*I+1.0D0)*X*PLM(I, I) |
| 1259 |
|
|
enddo |
| 1260 |
|
|
|
| 1261 |
|
|
DO I = 0, l |
| 1262 |
|
|
DO J = I+2, m |
| 1263 |
|
|
PLM(I, J)=((2.0D0*J-1.0D0)*X*PLM(I,J-1) - & |
| 1264 |
|
|
(I+J-1.0D0)*PLM(I,J-2))/(J-I) |
| 1265 |
|
|
end DO |
| 1266 |
|
|
end DO |
| 1267 |
|
|
|
| 1268 |
|
|
DLM(0, 0)=0.0D0 |
| 1269 |
|
|
|
| 1270 |
|
|
DO J = 1, m |
| 1271 |
|
|
DLM(0, J)=LS*J*(PLM(0,J-1)-X*PLM(0,J))/XS |
| 1272 |
|
|
end DO |
| 1273 |
|
|
|
| 1274 |
|
|
DO I = 1, l |
| 1275 |
|
|
DO J = I, m |
| 1276 |
|
|
DLM(I,J) = LS*I*X*PLM(I, J)/XS + (J+I)*(J-I+1.0D0)/XQ*PLM(I-1, J) |
| 1277 |
|
|
end DO |
| 1278 |
|
|
end DO |
| 1279 |
|
|
|
| 1280 |
|
|
RETURN |
| 1281 |
|
|
END SUBROUTINE Associated_Legendre |
| 1282 |
|
|
|
| 1283 |
|
|
|
| 1284 |
|
|
subroutine Orthogonal_Polynomial(x, m, function_type, pl, dpl) |
| 1285 |
|
|
|
| 1286 |
|
|
! Purpose: Compute orthogonal polynomials: Tn(x) or Un(x), |
| 1287 |
|
|
! or Ln(x) or Hn(x), and their derivatives |
| 1288 |
|
|
! Input : function_type --- Function code |
| 1289 |
|
|
! =1 for Chebyshev polynomial Tn(x) |
| 1290 |
|
|
! =2 for Chebyshev polynomial Un(x) |
| 1291 |
|
|
! =3 for Laguerre polynomial Ln(x) |
| 1292 |
|
|
! =4 for Hermite polynomial Hn(x) |
| 1293 |
|
|
! n --- Order of orthogonal polynomials |
| 1294 |
|
|
! x --- Argument of orthogonal polynomials |
| 1295 |
|
|
! Output: PL(n) --- Tn(x) or Un(x) or Ln(x) or Hn(x) |
| 1296 |
|
|
! DPL(n)--- Tn'(x) or Un'(x) or Ln'(x) or Hn'(x) |
| 1297 |
|
|
! |
| 1298 |
|
|
! adapted from the routines in |
| 1299 |
|
|
! COMPUTATION OF SPECIAL FUNCTIONS by Shanjie Zhang and Jianming Jin |
| 1300 |
|
|
! ISBN 0-471-11963-6 |
| 1301 |
|
|
! |
| 1302 |
|
|
! The original Fortran77 codes can be found here: |
| 1303 |
|
|
! http://iris-lee3.ece.uiuc.edu/~jjin/routines/routines.html |
| 1304 |
|
|
|
| 1305 |
|
|
real(kind=8), intent(in) :: x |
| 1306 |
|
|
integer, intent(in):: m |
| 1307 |
|
|
integer, intent(in):: function_type |
| 1308 |
|
|
real(kind=8), dimension(0:m), intent(inout) :: pl, dpl |
| 1309 |
|
|
|
| 1310 |
|
|
real(kind=8) :: a, b, c, y0, y1, dy0, dy1, yn, dyn |
| 1311 |
|
|
integer :: k |
| 1312 |
|
|
|
| 1313 |
|
|
A = 2.0D0 |
| 1314 |
|
|
B = 0.0D0 |
| 1315 |
|
|
C = 1.0D0 |
| 1316 |
|
|
Y0 = 1.0D0 |
| 1317 |
|
|
Y1 = 2.0D0*X |
| 1318 |
|
|
DY0 = 0.0D0 |
| 1319 |
|
|
DY1 = 2.0D0 |
| 1320 |
|
|
PL(0) = 1.0D0 |
| 1321 |
|
|
PL(1) = 2.0D0*X |
| 1322 |
|
|
DPL(0) = 0.0D0 |
| 1323 |
|
|
DPL(1) = 2.0D0 |
| 1324 |
|
|
IF (function_type.EQ.CHEBYSHEV_TN) THEN |
| 1325 |
|
|
Y1 = X |
| 1326 |
|
|
DY1 = 1.0D0 |
| 1327 |
|
|
PL(1) = X |
| 1328 |
|
|
DPL(1) = 1.0D0 |
| 1329 |
|
|
ELSE IF (function_type.EQ.LAGUERRE) THEN |
| 1330 |
|
|
Y1 = 1.0D0-X |
| 1331 |
|
|
DY1 = -1.0D0 |
| 1332 |
|
|
PL(1) = 1.0D0-X |
| 1333 |
|
|
DPL(1) = -1.0D0 |
| 1334 |
|
|
ENDIF |
| 1335 |
|
|
DO K = 2, m |
| 1336 |
|
|
IF (function_type.EQ.LAGUERRE) THEN |
| 1337 |
|
|
A = -1.0D0/K |
| 1338 |
|
|
B = 2.0D0+A |
| 1339 |
|
|
C = 1.0D0+A |
| 1340 |
|
|
ELSE IF (function_type.EQ.HERMITE) THEN |
| 1341 |
|
|
C = 2.0D0*(K-1.0D0) |
| 1342 |
|
|
ENDIF |
| 1343 |
|
|
YN = (A*X+B)*Y1-C*Y0 |
| 1344 |
|
|
DYN = A*Y1+(A*X+B)*DY1-C*DY0 |
| 1345 |
|
|
PL(K) = YN |
| 1346 |
|
|
DPL(K) = DYN |
| 1347 |
|
|
Y0 = Y1 |
| 1348 |
|
|
Y1 = YN |
| 1349 |
|
|
DY0 = DY1 |
| 1350 |
|
|
DY1 = DYN |
| 1351 |
|
|
end DO |
| 1352 |
|
|
RETURN |
| 1353 |
|
|
|
| 1354 |
|
|
end subroutine Orthogonal_Polynomial |
| 1355 |
|
|
|
| 1356 |
|
|
end module shapes |