ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/shapes.F90
Revision: 1717
Committed: Fri Nov 5 21:04:33 2004 UTC (19 years, 8 months ago) by chrisfen
File size: 50277 byte(s)
Log Message:
current status of the debugging process

File Contents

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