ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/shapes.F90
Revision: 1698
Committed: Tue Nov 2 15:28:25 2004 UTC (19 years, 8 months ago) by chrisfen
File size: 48133 byte(s)
Log Message:
Shapes looks like it's working

File Contents

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