ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_sticky_pair.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/calc_sticky_pair.F90 (file contents):
Revision 634 by gezelter, Tue Jul 15 17:10:50 2003 UTC vs.
Revision 635 by gezelter, Thu Jul 17 20:32:24 2003 UTC

# Line 9 | Line 9
9   !! @author Matthew Meineke
10   !! @author Christopher Fennel
11   !! @author J. Daniel Gezelter
12 < !! @version $Id: calc_sticky_pair.F90,v 1.10 2003-07-15 17:10:50 gezelter Exp $, $Date: 2003-07-15 17:10:50 $, $Name: not supported by cvs2svn $, $Revision: 1.10 $
12 > !! @version $Id: calc_sticky_pair.F90,v 1.11 2003-07-17 20:32:24 gezelter Exp $, $Date: 2003-07-17 20:32:24 $, $Name: not supported by cvs2svn $, $Revision: 1.11 $
13  
14   module sticky_pair
15  
# Line 27 | Line 27 | module sticky_pair
27    logical, save :: sticky_initialized = .false.
28    real( kind = dp ), save :: SSD_w0 = 0.0_dp
29    real( kind = dp ), save :: SSD_v0 = 0.0_dp
30 +  real( kind = dp ), save :: SSD_v0p = 0.0_dp
31    real( kind = dp ), save :: SSD_rl = 0.0_dp
32    real( kind = dp ), save :: SSD_ru = 0.0_dp
33 +  real( kind = dp ), save :: SSD_rlp = 0.0_dp
34    real( kind = dp ), save :: SSD_rup = 0.0_dp
35  
36    public :: check_sticky_FF
# Line 44 | Line 46 | contains
46      return
47    end subroutine check_sticky_FF
48  
49 <  subroutine set_sticky_params(sticky_w0, sticky_v0)
50 <    real( kind = dp ), intent(in) :: sticky_w0, sticky_v0
49 >  subroutine set_sticky_params(sticky_w0, sticky_v0, sticky_v0p, &
50 >       sticky_rl, sticky_ru, sticky_rlp, sticky_rup)
51 >
52 >    real( kind = dp ), intent(in) :: sticky_w0, sticky_v0, sticky_v0p
53 >    real( kind = dp ), intent(in) :: sticky_rl, sticky_ru
54 >    real( kind = dp ), intent(in) :: sticky_rlp, sticky_rup
55      
56      ! we could pass all 5 parameters if we felt like it...
57      
58      SSD_w0 = sticky_w0
59      SSD_v0 = sticky_v0
60 <    SSD_rl = 2.75_DP
61 <    SSD_ru = 3.35_DP
62 <    SSD_rup = 4.0_DP
60 >    SSD_v0p = sticky_v0p
61 >    SSD_rl = sticky_rl
62 >    SSD_ru = sticky_ru
63 >    SSD_rlp = sticky_rlp
64 >    SSD_rup = sticky_rup
65    
66      sticky_initialized = .true.
67      return
# Line 166 | Line 174 | contains
174  
175         if (do_pot) then
176   #ifdef IS_MPI
177 <          pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp)
178 <          pot_col(atom2) = pot_col(atom2) + 0.25d0*SSD_v0*(s*w + sp*wp)
177 >          pot_row(atom1) = pot_row(atom1) + 0.25d0*(SSD_v0*s*w + SSD_v0p*sp*wp)
178 >          pot_col(atom2) = pot_col(atom2) + 0.25d0*(SSD_v0*s*w + SSD_v0p*sp*wp)
179   #else
180 <          pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
180 >          pot = pot + 0.5d0*(SSD_v0*s*w + SSD_v0p*sp*wp)
181   #endif  
182         endif
183  
# Line 211 | Line 219 | contains
219         ! do the torques first since they are easy:
220         ! remember that these are still in the body fixed axes
221  
222 <       txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
223 <       tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
224 <       tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
222 >       txi = 0.5d0*(SSD_v0*s*dwidux + SSD_v0p*sp*dwipdux)
223 >       tyi = 0.5d0*(SSD_v0*s*dwiduy + SSD_v0p*sp*dwipduy)
224 >       tzi = 0.5d0*(SSD_v0*s*dwiduz + SSD_v0p*sp*dwipduz)
225  
226 <       txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
227 <       tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
228 <       tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
226 >       txj = 0.5d0*(SSD_v0*s*dwjdux + SSD_v0p*sp*dwjpdux)
227 >       tyj = 0.5d0*(SSD_v0*s*dwjduy + SSD_v0p*sp*dwjpduy)
228 >       tzj = 0.5d0*(SSD_v0*s*dwjduz + SSD_v0p*sp*dwjpduz)
229  
230         ! go back to lab frame using transpose of rotation matrix:
231  
# Line 248 | Line 256 | contains
256  
257         ! first rotate the i terms back into the lab frame:
258  
259 <       radcomxi = s*dwidx+sp*dwipdx
260 <       radcomyi = s*dwidy+sp*dwipdy
261 <       radcomzi = s*dwidz+sp*dwipdz
259 >       radcomxi = SSD_v0*s*dwidx+SSD_v0p*sp*dwipdx
260 >       radcomyi = SSD_v0*s*dwidy+SSD_v0p*sp*dwipdy
261 >       radcomzi = SSD_v0*s*dwidz+SSD_v0p*sp*dwipdz
262  
263 <       radcomxj = s*dwjdx+sp*dwjpdx
264 <       radcomyj = s*dwjdy+sp*dwjpdy
265 <       radcomzj = s*dwjdz+sp*dwjpdz
263 >       radcomxj = SSD_v0*s*dwjdx+SSD_v0p*sp*dwjpdx
264 >       radcomyj = SSD_v0*s*dwjdy+SSD_v0p*sp*dwjpdy
265 >       radcomzj = SSD_v0*s*dwjdz+SSD_v0p*sp*dwjpdz
266  
267   #ifdef IS_MPI    
268         fxii = a_Row(1,atom1)*(radcomxi) + &
# Line 308 | Line 316 | contains
316  
317         ! now assemble these with the radial-only terms:
318  
319 <       fxradial = 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + fxii + fxji)
320 <       fyradial = 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + fyii + fyji)
321 <       fzradial = 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + fzii + fzji)
319 >       fxradial = 0.5d0*(SSD_v0*dsdr*drdx*w + SSD_v0p*dspdr*drdx*wp + fxii + fxji)
320 >       fyradial = 0.5d0*(SSD_v0*dsdr*drdy*w + SSD_v0p*dspdr*drdy*wp + fyii + fyji)
321 >       fzradial = 0.5d0*(SSD_v0*dsdr*drdz*w + SSD_v0p*dspdr*drdz*wp + fzii + fzji)
322  
323   #ifdef IS_MPI
324         f_Row(1,atom1) = f_Row(1,atom1) + fxradial
# Line 357 | Line 365 | contains
365  
366    !! calculates the switching functions and their derivatives for a given
367    subroutine calc_sw_fnc(r, s, sp, dsdr, dspdr)
368 <          
368 >    
369      real (kind=dp), intent(in) :: r
370      real (kind=dp), intent(inout) :: s, sp, dsdr, dspdr
371 <
371 >    
372      ! distances must be in angstroms
373      
374      if (r.lt.SSD_rl) then
375         s = 1.0d0
368       sp = 1.0d0
376         dsdr = 0.0d0
377 +    elseif (r.gt.SSD_ru) then
378 +       s = 0.0d0
379 +       dsdr = 0.0d0
380 +    else
381 +       s = ((SSD_ru + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_ru-r)**2) / &
382 +            ((SSD_ru - SSD_rl)**3)
383 +       dsdr = 6.0d0*(r-SSD_ru)*(r-SSD_rl)/((SSD_ru - SSD_rl)**3)
384 +    endif
385 +
386 +    if (r.lt.SSD_rlp) then
387 +       sp = 1.0d0      
388         dspdr = 0.0d0
389      elseif (r.gt.SSD_rup) then
372       s = 0.0d0
390         sp = 0.0d0
374       dsdr = 0.0d0
391         dspdr = 0.0d0
392      else
393 <       sp = ((SSD_rup + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_rup-r)**2) / &
394 <            ((SSD_rup - SSD_rl)**3)
395 <       dspdr = 6.0d0*(r-SSD_rup)*(r-SSD_rl)/((SSD_rup - SSD_rl)**3)
380 <      
381 <       if (r.gt.SSD_ru) then
382 <          s = 0.0d0
383 <          dsdr = 0.0d0
384 <       else
385 <          s = ((SSD_ru + 2.0d0*r - 3.0d0*SSD_rl) * (SSD_ru-r)**2) / &
386 <               ((SSD_ru - SSD_rl)**3)
387 <          dsdr = 6.0d0*(r-SSD_ru)*(r-SSD_rl)/((SSD_ru - SSD_rl)**3)
388 <       endif
393 >       sp = ((SSD_rup + 2.0d0*r - 3.0d0*SSD_rlp) * (SSD_rup-r)**2) / &
394 >            ((SSD_rup - SSD_rlp)**3)
395 >       dspdr = 6.0d0*(r-SSD_rup)*(r-SSD_rlp)/((SSD_rup - SSD_rlp)**3)      
396      endif
397      
398      return

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines