ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/shapes.F90
Revision: 2214
Committed: Wed Apr 27 20:14:03 2005 UTC (19 years, 2 months ago) by chrisfen
File size: 53824 byte(s)
Log Message:
Got rid of write statements and am closer to a working shapes

File Contents

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