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