ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/shapes.F90
Revision: 1948
Committed: Fri Jan 14 20:31:16 2005 UTC (19 years, 5 months ago) by gezelter
File size: 50341 byte(s)
Log Message:
separating modules and C/Fortran interface subroutines

File Contents

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