ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/shapes.F90
Revision: 1608
Committed: Wed Oct 20 04:02:48 2004 UTC (19 years, 9 months ago) by gezelter
File size: 49564 byte(s)
Log Message:
name sanity on the fortran side

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