ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/shapes.F90
Revision: 1464
Committed: Fri Jul 9 19:29:05 2010 UTC (15 years, 1 month ago) by gezelter
File size: 50882 byte(s)
Log Message:
removing cruft (atom numbers, do_pot, do_stress) from many modules and
force managers

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

Properties

Name Value
svn:keywords Author Id Revision Date