ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/shapes.F90
Revision: 1707
Committed: Thu Nov 4 16:20:37 2004 UTC (19 years, 8 months ago) by gezelter
File size: 50512 byte(s)
Log Message:
Breaky Breaky

File Contents

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