ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/shapes.F90
Revision: 1386
Committed: Fri Oct 23 18:41:09 2009 UTC (15 years, 9 months ago) by gezelter
File size: 51044 byte(s)
Log Message:
removing MPI responsibilities from the lowest level force routines.  This is
in preparation for migrating these routines (LJ, electrostatic, eam, suttonchen,
gay-berne, sticky) to C++

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