ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90
(Generate patch)

Comparing trunk/OOPSE-4/src/UseTheForce/DarkSide/electrostatic.F90 (file contents):
Revision 2296 by chrisfen, Thu Sep 15 00:13:56 2005 UTC vs.
Revision 2302 by chrisfen, Fri Sep 16 16:07:39 2005 UTC

# Line 54 | Line 54 | module electrostatic_module
54  
55    PRIVATE
56  
57 + #define __FORTRAN90
58 + #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 +
60    !! these prefactors convert the multipole interactions into kcal / mol
61    !! all were computed assuming distances are measured in angstroms
62    !! Charge-Charge, assuming charges are measured in electrons
# Line 68 | Line 71 | module electrostatic_module
71    !! This unit is also known affectionately as an esu centi-barn.
72    real(kind=dp), parameter :: pre14 = 69.13373_dp
73  
74 +  !! variables to handle different summation methods for long-range electrostatics:
75 +  integer, save :: summationMethod = NONE
76 +  logical, save :: summationMethodChecked = .false.
77 +  real(kind=DP), save :: defaultCutoff = 0.0_DP
78 +  logical, save :: haveDefaultCutoff = .false.
79 +  real(kind=DP), save :: dampingAlpha = 0.0_DP
80 +  logical, save :: haveDampingAlpha = .false.
81 +  real(kind=DP), save :: dielectric = 0.0_DP
82 +  logical, save :: haveDielectric = .false.
83 +  real(kind=DP), save :: constERFC = 0.0_DP
84 +  real(kind=DP), save :: constEXP = 0.0_DP
85 +  logical, save :: haveDWAconstants = .false.
86 +
87 +
88 +  public :: setElectrostaticSummationMethod
89 +  public :: setElectrostaticCutoffRadius
90 +  public :: setDampedWolfAlpha
91 +  public :: setReactionFieldDielectric
92    public :: newElectrostaticType
93    public :: setCharge
94    public :: setDipoleMoment
# Line 96 | Line 117 | contains
117  
118   contains
119  
120 +  subroutine setElectrostaticSummationMethod(the_ESM)
121 +
122 +    integer, intent(in) :: the_ESM    
123 +
124 +    if ((the_ESM .le. 0) .or. (the_ESM .gt. REACTION_FIELD)) then
125 +       call handleError("setElectrostaticSummationMethod", "Unsupported Summation Method")
126 +    endif
127 +
128 +  end subroutine setElectrostaticSummationMethod
129 +
130 +  subroutine setElectrostaticCutoffRadius(thisRcut)
131 +    real(kind=dp), intent(in) :: thisRcut
132 +    defaultCutoff = thisRcut
133 +    haveDefaultCutoff = .true.
134 +  end subroutine setElectrostaticCutoffRadius
135 +
136 +  subroutine setDampedWolfAlpha(thisAlpha)
137 +    real(kind=dp), intent(in) :: thisAlpha
138 +    dampingAlpha = thisAlpha
139 +    haveDampingAlpha = .true.
140 +  end subroutine setDampedWolfAlpha
141 +  
142 +  subroutine setReactionFieldDielectric(thisDielectric)
143 +    real(kind=dp), intent(in) :: thisDielectric
144 +    dielectric = thisDielectric
145 +    haveDielectric = .true.
146 +  end subroutine setReactionFieldDielectric
147 +
148    subroutine newElectrostaticType(c_ident, is_Charge, is_Dipole, &
149         is_SplitDipole, is_Quadrupole, is_Tap, status)
150  
# Line 305 | Line 354 | contains
354  
355      dm = ElectrostaticMap(atid)%dipole_moment
356    end function getDipoleMoment
357 +
358 +  subroutine checkSummationMethod()
359 +
360 +    if (summationMethod .eq. DAMPED_WOLF) then
361 +       if (.not.haveDWAconstants) then
362 +          
363 +          if (.not.haveDampingAlpha) then
364 +             call handleError("checkSummationMethod", "no Damping Alpha set!")
365 +          endif
366 +          
367 +          if (.not.haveDefaultCutoff) then
368 +             call handleError("checkSummationMethod", "no Default Cutoff set!")
369 +          endif
370 +
371 +          constEXP = exp(-dampingAlpha*dampingAlpha*defaultCutoff*defaultCutoff)
372 +          constERFC = erfc(dampingAlpha*defaultCutoff)
373 +          
374 +          haveDWAconstants = .true.
375 +       endif
376 +    endif
377 +
378 +    if (summationMethod .eq. REACTION_FIELD) then
379 +       if (.not.haveDielectric) then
380 +          call handleError("checkSummationMethod", "no reaction field Dielectric set!")
381 +       endif
382 +    endif
383 +
384 +    summationMethodChecked = .true.
385 +  end subroutine checkSummationMethod
386  
387 +
388 +
389    subroutine doElectrostaticPair(atom1, atom2, d, rij, r2, sw, &
390         vpair, fpair, pot, eFrame, f, t, do_pot, corrMethod, rcuti)
391  
# Line 353 | Line 433 | contains
433         return
434      end if
435  
436 +    if (.not.summationMethodChecked) then
437 +       call checkSummationMethod()
438 +    endif
439 +
440 +
441   #ifdef IS_MPI
442      me1 = atid_Row(atom1)
443      me2 = atid_Col(atom2)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines