ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/sticky.F90
Revision: 2189
Committed: Wed Apr 13 20:36:45 2005 UTC (19 years, 3 months ago) by chuckv
File size: 17397 byte(s)
Log Message:
Added destroy methods for Fortran modules.

File Contents

# Content
1 !!
2 !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 !!
4 !! The University of Notre Dame grants you ("Licensee") a
5 !! non-exclusive, royalty free, license to use, modify and
6 !! redistribute this software in source and binary code form, provided
7 !! that the following conditions are met:
8 !!
9 !! 1. Acknowledgement of the program authors must be made in any
10 !! publication of scientific results based in part on use of the
11 !! program. An acceptable form of acknowledgement is citation of
12 !! the article in which the program was described (Matthew
13 !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 !! Parallel Simulation Engine for Molecular Dynamics,"
16 !! J. Comput. Chem. 26, pp. 252-271 (2005))
17 !!
18 !! 2. Redistributions of source code must retain the above copyright
19 !! notice, this list of conditions and the following disclaimer.
20 !!
21 !! 3. Redistributions in binary form must reproduce the above copyright
22 !! notice, this list of conditions and the following disclaimer in the
23 !! documentation and/or other materials provided with the
24 !! distribution.
25 !!
26 !! This software is provided "AS IS," without a warranty of any
27 !! kind. All express or implied conditions, representations and
28 !! warranties, including any implied warranty of merchantability,
29 !! fitness for a particular purpose or non-infringement, are hereby
30 !! excluded. The University of Notre Dame and its licensors shall not
31 !! be liable for any damages suffered by licensee as a result of
32 !! using, modifying or distributing the software or its
33 !! derivatives. In no event will the University of Notre Dame or its
34 !! licensors be liable for any lost revenue, profit or data, or for
35 !! direct, indirect, special, consequential, incidental or punitive
36 !! damages, however caused and regardless of the theory of liability,
37 !! arising out of the use of or inability to use software, even if the
38 !! University of Notre Dame has been advised of the possibility of
39 !! such damages.
40 !!
41
42 !! This Module Calculates forces due to SSD potential and VDW interactions
43 !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
44
45 !! This module contains the Public procedures:
46
47
48 !! Corresponds to the force field defined in ssd_FF.cpp
49 !! @author Charles F. Vardeman II
50 !! @author Matthew Meineke
51 !! @author Christopher Fennell
52 !! @author J. Daniel Gezelter
53 !! @version $Id: sticky.F90,v 1.6 2005-04-13 20:36:45 chuckv Exp $, $Date: 2005-04-13 20:36:45 $, $Name: not supported by cvs2svn $, $Revision: 1.6 $
54
55 module sticky
56
57 use force_globals
58 use definitions
59 use atype_module
60 use vector_class
61 use simulation
62 use status
63 #ifdef IS_MPI
64 use mpiSimulation
65 #endif
66 implicit none
67
68 PRIVATE
69
70 public :: newStickyType
71 public :: do_sticky_pair
72 public :: destroyStickyTypes
73
74
75 type :: StickyList
76 integer :: c_ident
77 real( kind = dp ) :: w0 = 0.0_dp
78 real( kind = dp ) :: v0 = 0.0_dp
79 real( kind = dp ) :: v0p = 0.0_dp
80 real( kind = dp ) :: rl = 0.0_dp
81 real( kind = dp ) :: ru = 0.0_dp
82 real( kind = dp ) :: rlp = 0.0_dp
83 real( kind = dp ) :: rup = 0.0_dp
84 real( kind = dp ) :: rbig = 0.0_dp
85 end type StickyList
86
87 type(StickyList), dimension(:),allocatable :: StickyMap
88
89 contains
90
91 subroutine newStickyType(c_ident, w0, v0, v0p, rl, ru, rlp, rup, isError)
92
93 integer, intent(in) :: c_ident
94 integer, intent(inout) :: isError
95 real( kind = dp ), intent(in) :: w0, v0, v0p
96 real( kind = dp ), intent(in) :: rl, ru
97 real( kind = dp ), intent(in) :: rlp, rup
98 integer :: nATypes, myATID
99
100
101 isError = 0
102 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
103
104 !! Be simple-minded and assume that we need a StickyMap that
105 !! is the same size as the total number of atom types
106
107 if (.not.allocated(StickyMap)) then
108
109 nAtypes = getSize(atypes)
110
111 if (nAtypes == 0) then
112 isError = -1
113 return
114 end if
115
116 if (.not. allocated(StickyMap)) then
117 allocate(StickyMap(nAtypes))
118 endif
119
120 end if
121
122 if (myATID .gt. size(StickyMap)) then
123 isError = -1
124 return
125 endif
126
127 ! set the values for StickyMap for this atom type:
128
129 StickyMap(myATID)%c_ident = c_ident
130
131 ! we could pass all 5 parameters if we felt like it...
132
133 StickyMap(myATID)%w0 = w0
134 StickyMap(myATID)%v0 = v0
135 StickyMap(myATID)%v0p = v0p
136 StickyMap(myATID)%rl = rl
137 StickyMap(myATID)%ru = ru
138 StickyMap(myATID)%rlp = rlp
139 StickyMap(myATID)%rup = rup
140
141 if (StickyMap(myATID)%ru .gt. StickyMap(myATID)%rup) then
142 StickyMap(myATID)%rbig = StickyMap(myATID)%ru
143 else
144 StickyMap(myATID)%rbig = StickyMap(myATID)%rup
145 endif
146
147 return
148 end subroutine newStickyType
149
150 subroutine do_sticky_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
151 pot, A, f, t, do_pot)
152
153 !! This routine does only the sticky portion of the SSD potential
154 !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
155 !! The Lennard-Jones and dipolar interaction must be handled separately.
156
157 !! We assume that the rotation matrices have already been calculated
158 !! and placed in the A array.
159
160 !! i and j are pointers to the two SSD atoms
161
162 integer, intent(in) :: atom1, atom2
163 real (kind=dp), intent(inout) :: rij, r2
164 real (kind=dp), dimension(3), intent(in) :: d
165 real (kind=dp), dimension(3), intent(inout) :: fpair
166 real (kind=dp) :: pot, vpair, sw
167 real (kind=dp), dimension(9,nLocal) :: A
168 real (kind=dp), dimension(3,nLocal) :: f
169 real (kind=dp), dimension(3,nLocal) :: t
170 logical, intent(in) :: do_pot
171
172 real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
173 real (kind=dp) :: r3, r5, r6, s, sp, dsdr, dspdr
174 real (kind=dp) :: wi, wj, w, wip, wjp, wp
175 real (kind=dp) :: dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
176 real (kind=dp) :: dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
177 real (kind=dp) :: dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
178 real (kind=dp) :: dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
179 real (kind=dp) :: zif, zis, zjf, zjs, uglyi, uglyj
180 real (kind=dp) :: drdx, drdy, drdz
181 real (kind=dp) :: txi, tyi, tzi, txj, tyj, tzj
182 real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj
183 real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji
184 real (kind=dp) :: fxradial, fyradial, fzradial
185 real (kind=dp) :: rijtest, rjitest
186 real (kind=dp) :: radcomxi, radcomyi, radcomzi
187 real (kind=dp) :: radcomxj, radcomyj, radcomzj
188 integer :: id1, id2
189 integer :: me1, me2
190 real (kind=dp) :: w0, v0, v0p, rl, ru, rlp, rup, rbig
191
192 if (.not.allocated(StickyMap)) then
193 call handleError("sticky", "no StickyMap was present before first call of do_sticky_pair!")
194 return
195 end if
196
197 #ifdef IS_MPI
198 me1 = atid_Row(atom1)
199 me2 = atid_Col(atom2)
200 #else
201 me1 = atid(atom1)
202 me2 = atid(atom2)
203 #endif
204
205 if (me1.eq.me2) then
206 w0 = StickyMap(me1)%w0
207 v0 = StickyMap(me1)%v0
208 v0p = StickyMap(me1)%v0p
209 rl = StickyMap(me1)%rl
210 ru = StickyMap(me1)%ru
211 rlp = StickyMap(me1)%rlp
212 rup = StickyMap(me1)%rup
213 rbig = StickyMap(me1)%rbig
214 else
215 ! This is silly, but if you want 2 sticky types in your
216 ! simulation, we'll let you do it with the Lorentz-
217 ! Berthelot mixing rules.
218 ! (Warning: you'll be SLLLLLLLLLLLLLLLOOOOOOOOOOWWWWWWWWWWW)
219 rl = 0.5_dp * ( StickyMap(me1)%rl + StickyMap(me2)%rl )
220 ru = 0.5_dp * ( StickyMap(me1)%ru + StickyMap(me2)%ru )
221 rlp = 0.5_dp * ( StickyMap(me1)%rlp + StickyMap(me2)%rlp )
222 rup = 0.5_dp * ( StickyMap(me1)%rup + StickyMap(me2)%rup )
223 rbig = max(ru, rup)
224 w0 = sqrt( StickyMap(me1)%w0 * StickyMap(me2)%w0 )
225 v0 = sqrt( StickyMap(me1)%v0 * StickyMap(me2)%v0 )
226 v0p = sqrt( StickyMap(me1)%v0p * StickyMap(me2)%v0p )
227 endif
228
229 if ( rij .LE. rbig ) then
230
231 r3 = r2*rij
232 r5 = r3*r2
233
234 drdx = d(1) / rij
235 drdy = d(2) / rij
236 drdz = d(3) / rij
237
238 #ifdef IS_MPI
239 ! rotate the inter-particle separation into the two different
240 ! body-fixed coordinate systems:
241
242 xi = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
243 yi = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
244 zi = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
245
246 ! negative sign because this is the vector from j to i:
247
248 xj = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
249 yj = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
250 zj = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
251 #else
252 ! rotate the inter-particle separation into the two different
253 ! body-fixed coordinate systems:
254
255 xi = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
256 yi = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
257 zi = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
258
259 ! negative sign because this is the vector from j to i:
260
261 xj = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
262 yj = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
263 zj = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
264 #endif
265
266 xi2 = xi*xi
267 yi2 = yi*yi
268 zi2 = zi*zi
269
270 xj2 = xj*xj
271 yj2 = yj*yj
272 zj2 = zj*zj
273
274 call calc_sw_fnc(rij, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
275
276 wi = 2.0d0*(xi2-yi2)*zi / r3
277 wj = 2.0d0*(xj2-yj2)*zj / r3
278 w = wi+wj
279
280 zif = zi/rij - 0.6d0
281 zis = zi/rij + 0.8d0
282
283 zjf = zj/rij - 0.6d0
284 zjs = zj/rij + 0.8d0
285
286 wip = zif*zif*zis*zis - w0
287 wjp = zjf*zjf*zjs*zjs - w0
288 wp = wip + wjp
289
290 vpair = vpair + 0.5d0*(v0*s*w + v0p*sp*wp)
291 if (do_pot) then
292 #ifdef IS_MPI
293 pot_row(atom1) = pot_row(atom1) + 0.25d0*(v0*s*w + v0p*sp*wp)*sw
294 pot_col(atom2) = pot_col(atom2) + 0.25d0*(v0*s*w + v0p*sp*wp)*sw
295 #else
296 pot = pot + 0.5d0*(v0*s*w + v0p*sp*wp)*sw
297 #endif
298 endif
299
300 dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5
301 dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5
302 dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5
303
304 dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5
305 dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5
306 dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5
307
308 uglyi = zif*zif*zis + zif*zis*zis
309 uglyj = zjf*zjf*zjs + zjf*zjs*zjs
310
311 dwipdx = -2.0d0*xi*zi*uglyi/r3
312 dwipdy = -2.0d0*yi*zi*uglyi/r3
313 dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
314
315 dwjpdx = -2.0d0*xj*zj*uglyj/r3
316 dwjpdy = -2.0d0*yj*zj*uglyj/r3
317 dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
318
319 dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
320 dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
321 dwiduz = - 8.0d0*xi*yi*zi/r3
322
323 dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
324 dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
325 dwjduz = - 8.0d0*xj*yj*zj/r3
326
327 dwipdux = 2.0d0*yi*uglyi/rij
328 dwipduy = -2.0d0*xi*uglyi/rij
329 dwipduz = 0.0d0
330
331 dwjpdux = 2.0d0*yj*uglyj/rij
332 dwjpduy = -2.0d0*xj*uglyj/rij
333 dwjpduz = 0.0d0
334
335 ! do the torques first since they are easy:
336 ! remember that these are still in the body fixed axes
337
338 txi = 0.5d0*(v0*s*dwidux + v0p*sp*dwipdux)*sw
339 tyi = 0.5d0*(v0*s*dwiduy + v0p*sp*dwipduy)*sw
340 tzi = 0.5d0*(v0*s*dwiduz + v0p*sp*dwipduz)*sw
341
342 txj = 0.5d0*(v0*s*dwjdux + v0p*sp*dwjpdux)*sw
343 tyj = 0.5d0*(v0*s*dwjduy + v0p*sp*dwjpduy)*sw
344 tzj = 0.5d0*(v0*s*dwjduz + v0p*sp*dwjpduz)*sw
345
346 ! go back to lab frame using transpose of rotation matrix:
347
348 #ifdef IS_MPI
349 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*txi + &
350 a_Row(4,atom1)*tyi + a_Row(7,atom1)*tzi
351 t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*txi + &
352 a_Row(5,atom1)*tyi + a_Row(8,atom1)*tzi
353 t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*txi + &
354 a_Row(6,atom1)*tyi + a_Row(9,atom1)*tzi
355
356 t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*txj + &
357 a_Col(4,atom2)*tyj + a_Col(7,atom2)*tzj
358 t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*txj + &
359 a_Col(5,atom2)*tyj + a_Col(8,atom2)*tzj
360 t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*txj + &
361 a_Col(6,atom2)*tyj + a_Col(9,atom2)*tzj
362 #else
363 t(1,atom1) = t(1,atom1) + a(1,atom1)*txi + a(4,atom1)*tyi + a(7,atom1)*tzi
364 t(2,atom1) = t(2,atom1) + a(2,atom1)*txi + a(5,atom1)*tyi + a(8,atom1)*tzi
365 t(3,atom1) = t(3,atom1) + a(3,atom1)*txi + a(6,atom1)*tyi + a(9,atom1)*tzi
366
367 t(1,atom2) = t(1,atom2) + a(1,atom2)*txj + a(4,atom2)*tyj + a(7,atom2)*tzj
368 t(2,atom2) = t(2,atom2) + a(2,atom2)*txj + a(5,atom2)*tyj + a(8,atom2)*tzj
369 t(3,atom2) = t(3,atom2) + a(3,atom2)*txj + a(6,atom2)*tyj + a(9,atom2)*tzj
370 #endif
371 ! Now, on to the forces:
372
373 ! first rotate the i terms back into the lab frame:
374
375 radcomxi = (v0*s*dwidx+v0p*sp*dwipdx)*sw
376 radcomyi = (v0*s*dwidy+v0p*sp*dwipdy)*sw
377 radcomzi = (v0*s*dwidz+v0p*sp*dwipdz)*sw
378
379 radcomxj = (v0*s*dwjdx+v0p*sp*dwjpdx)*sw
380 radcomyj = (v0*s*dwjdy+v0p*sp*dwjpdy)*sw
381 radcomzj = (v0*s*dwjdz+v0p*sp*dwjpdz)*sw
382
383 #ifdef IS_MPI
384 fxii = a_Row(1,atom1)*(radcomxi) + &
385 a_Row(4,atom1)*(radcomyi) + &
386 a_Row(7,atom1)*(radcomzi)
387 fyii = a_Row(2,atom1)*(radcomxi) + &
388 a_Row(5,atom1)*(radcomyi) + &
389 a_Row(8,atom1)*(radcomzi)
390 fzii = a_Row(3,atom1)*(radcomxi) + &
391 a_Row(6,atom1)*(radcomyi) + &
392 a_Row(9,atom1)*(radcomzi)
393
394 fxjj = a_Col(1,atom2)*(radcomxj) + &
395 a_Col(4,atom2)*(radcomyj) + &
396 a_Col(7,atom2)*(radcomzj)
397 fyjj = a_Col(2,atom2)*(radcomxj) + &
398 a_Col(5,atom2)*(radcomyj) + &
399 a_Col(8,atom2)*(radcomzj)
400 fzjj = a_Col(3,atom2)*(radcomxj)+ &
401 a_Col(6,atom2)*(radcomyj) + &
402 a_Col(9,atom2)*(radcomzj)
403 #else
404 fxii = a(1,atom1)*(radcomxi) + &
405 a(4,atom1)*(radcomyi) + &
406 a(7,atom1)*(radcomzi)
407 fyii = a(2,atom1)*(radcomxi) + &
408 a(5,atom1)*(radcomyi) + &
409 a(8,atom1)*(radcomzi)
410 fzii = a(3,atom1)*(radcomxi) + &
411 a(6,atom1)*(radcomyi) + &
412 a(9,atom1)*(radcomzi)
413
414 fxjj = a(1,atom2)*(radcomxj) + &
415 a(4,atom2)*(radcomyj) + &
416 a(7,atom2)*(radcomzj)
417 fyjj = a(2,atom2)*(radcomxj) + &
418 a(5,atom2)*(radcomyj) + &
419 a(8,atom2)*(radcomzj)
420 fzjj = a(3,atom2)*(radcomxj)+ &
421 a(6,atom2)*(radcomyj) + &
422 a(9,atom2)*(radcomzj)
423 #endif
424
425 fxij = -fxii
426 fyij = -fyii
427 fzij = -fzii
428
429 fxji = -fxjj
430 fyji = -fyjj
431 fzji = -fzjj
432
433 ! now assemble these with the radial-only terms:
434
435 fxradial = 0.5d0*(v0*dsdr*drdx*w + v0p*dspdr*drdx*wp + fxii + fxji)
436 fyradial = 0.5d0*(v0*dsdr*drdy*w + v0p*dspdr*drdy*wp + fyii + fyji)
437 fzradial = 0.5d0*(v0*dsdr*drdz*w + v0p*dspdr*drdz*wp + fzii + fzji)
438
439 #ifdef IS_MPI
440 f_Row(1,atom1) = f_Row(1,atom1) + fxradial
441 f_Row(2,atom1) = f_Row(2,atom1) + fyradial
442 f_Row(3,atom1) = f_Row(3,atom1) + fzradial
443
444 f_Col(1,atom2) = f_Col(1,atom2) - fxradial
445 f_Col(2,atom2) = f_Col(2,atom2) - fyradial
446 f_Col(3,atom2) = f_Col(3,atom2) - fzradial
447 #else
448 f(1,atom1) = f(1,atom1) + fxradial
449 f(2,atom1) = f(2,atom1) + fyradial
450 f(3,atom1) = f(3,atom1) + fzradial
451
452 f(1,atom2) = f(1,atom2) - fxradial
453 f(2,atom2) = f(2,atom2) - fyradial
454 f(3,atom2) = f(3,atom2) - fzradial
455 #endif
456
457 #ifdef IS_MPI
458 id1 = AtomRowToGlobal(atom1)
459 id2 = AtomColToGlobal(atom2)
460 #else
461 id1 = atom1
462 id2 = atom2
463 #endif
464
465 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
466
467 fpair(1) = fpair(1) + fxradial
468 fpair(2) = fpair(2) + fyradial
469 fpair(3) = fpair(3) + fzradial
470
471 endif
472 endif
473 end subroutine do_sticky_pair
474
475 !! calculates the switching functions and their derivatives for a given
476 subroutine calc_sw_fnc(r, rl, ru, rlp, rup, s, sp, dsdr, dspdr)
477
478 real (kind=dp), intent(in) :: r, rl, ru, rlp, rup
479 real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
480
481 ! distances must be in angstroms
482
483 if (r.lt.rl) then
484 s = 1.0d0
485 dsdr = 0.0d0
486 elseif (r.gt.ru) then
487 s = 0.0d0
488 dsdr = 0.0d0
489 else
490 s = ((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) / &
491 ((ru - rl)**3)
492 dsdr = 6.0d0*(r-ru)*(r-rl)/((ru - rl)**3)
493 endif
494
495 if (r.lt.rlp) then
496 sp = 1.0d0
497 dspdr = 0.0d0
498 elseif (r.gt.rup) then
499 sp = 0.0d0
500 dspdr = 0.0d0
501 else
502 sp = ((rup + 2.0d0*r - 3.0d0*rlp) * (rup-r)**2) / &
503 ((rup - rlp)**3)
504 dspdr = 6.0d0*(r-rup)*(r-rlp)/((rup - rlp)**3)
505 endif
506
507 return
508 end subroutine calc_sw_fnc
509
510 subroutine destroyStickyTypes()
511 if(allocated(StickyMap)) deallocate(StickyMap)
512 end subroutine destroyStickyTypes
513 end module sticky