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