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

Comparing trunk/OOPSE-3.0/src/UseTheForce/doForces.F90 (file contents):
Revision 2290 by gezelter, Wed Sep 14 20:28:05 2005 UTC vs.
Revision 2295 by chrisfen, Thu Sep 15 00:13:15 2005 UTC

# Line 45 | Line 45
45  
46   !! @author Charles F. Vardeman II
47   !! @author Matthew Meineke
48 < !! @version $Id: doForces.F90,v 1.41 2005-09-14 20:28:05 gezelter Exp $, $Date: 2005-09-14 20:28:05 $, $Name: not supported by cvs2svn $, $Revision: 1.41 $
48 > !! @version $Id: doForces.F90,v 1.42 2005-09-15 00:13:14 chrisfen Exp $, $Date: 2005-09-15 00:13:14 $, $Name: not supported by cvs2svn $, $Revision: 1.42 $
49  
50  
51   module doForces
# Line 58 | Line 58 | module doForces
58    use lj
59    use sticky
60    use electrostatic_module
61 <  use reaction_field
61 >  use reaction_field_module
62    use gb_pair
63    use shapes
64    use vector_class
# Line 74 | Line 74 | module doForces
74   #define __FORTRAN90
75   #include "UseTheForce/fSwitchingFunction.h"
76   #include "UseTheForce/fCutoffPolicy.h"
77 + #include "UseTheForce/fCoulombicCorrection.h"
78   #include "UseTheForce/DarkSide/fInteractionMap.h"
79  
80  
# Line 135 | Line 136 | module doForces
136    type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap
137  
138    integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY
139 +  integer, save :: coulombicCorrection = NONE
140    real(kind=dp),save :: defaultRcut, defaultRsw, defaultRlist
141 +  real(kind=dp),save :: rcuti
142    
143   contains
144  
# Line 408 | Line 411 | contains
411      enddo
412      
413      haveGtypeCutoffMap = .true.
414 <    
412 <  end subroutine createGtypeCutoffMap
413 <  
414 <  subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
415 <    real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
416 <    integer, intent(in) :: cutPolicy
417 <    
418 <    defaultRcut = defRcut
419 <    defaultRsw = defRsw
420 <    defaultRlist = defRlist
421 <    cutoffPolicy = cutPolicy
422 <  end subroutine setDefaultCutoffs
423 <  
424 <  subroutine setCutoffPolicy(cutPolicy)
414 >   end subroutine createGtypeCutoffMap
415  
416 +   subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
417 +     real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
418       integer, intent(in) :: cutPolicy
419 +
420 +     defaultRcut = defRcut
421 +     defaultRsw = defRsw
422 +     defaultRlist = defRlist
423       cutoffPolicy = cutPolicy
424 <     call createGtypeCutoffMap()
424 >     rcuti = 1.0_dp / defaultRcut
425 >   end subroutine setDefaultCutoffs
426  
427 +   subroutine setCutoffPolicy(cutPolicy)
428 +
429 +     integer, intent(in) :: cutPolicy
430 +     cutoffPolicy = cutPolicy
431 +     call createGtypeCutoffMap()
432     end subroutine setCutoffPolicy
433      
434      
435    subroutine setSimVariables()
436      SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
437      SIM_uses_EAM = SimUsesEAM()
436    SIM_uses_RF = SimUsesRF()
438      SIM_requires_postpair_calc = SimRequiresPostpairCalc()
439      SIM_requires_prepair_calc = SimRequiresPrepairCalc()
440      SIM_uses_PBC = SimUsesPBC()
441 +    SIM_uses_RF = SimUsesRF()
442  
443      haveSIMvariables = .true.
444  
# Line 503 | Line 505 | contains
505    end subroutine doReadyCheck
506  
507  
508 <  subroutine init_FF(use_RF, use_UW, use_DW, thisStat)
508 >  subroutine init_FF(use_RF, correctionMethod, dampingAlpha, thisStat)
509  
510      logical, intent(in) :: use_RF
511 <    logical, intent(in) :: use_UW
512 <    logical, intent(in) :: use_DW
511 >    integer, intent(in) :: correctionMethod
512 >    real(kind=dp), intent(in) :: dampingAlpha
513      integer, intent(out) :: thisStat  
514      integer :: my_status, nMatches
513    integer :: corrMethod
515      integer, pointer :: MatchList(:) => null()
516      real(kind=dp) :: rcut, rrf, rt, dielect
517  
# Line 521 | Line 522 | contains
522      FF_uses_RF = use_RF
523  
524      !! set the electrostatic correction method
525 <    if (use_UW) then
525 >    select case(coulombicCorrection)
526 >    case(NONE)
527 >       corrMethod = 0
528 >    case(UNDAMPED_WOLF)
529         corrMethod = 1
530 <    elseif (use_DW) then
530 >    case(WOLF)
531         corrMethod = 2
532 <    else
533 <       corrMethod = 0
534 <    endif
535 <    
532 >    case (REACTION_FIELD)
533 >       corrMethod = 3
534 >    case default
535 >       call handleError("init_FF", "Unknown Coulombic Correction Method")
536 >       return
537 >    end select
538 >        
539      !! init_FF is called *after* all of the atom types have been
540      !! defined in atype_module using the new_atype subroutine.
541      !!
# Line 566 | Line 573 | contains
573            call initialize_rf(dielect)
574         endif
575      else
576 <       if (FF_uses_RF) then          
576 >       if ((corrMethod == 3) .or. FF_uses_RF) then
577            write(default_error,*) 'Using Reaction Field with no dipoles?  Huh?'
578            thisStat = -1
579            haveSaneForceField = .false.
# Line 984 | Line 991 | contains
991  
992      if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
993  
994 <       if (FF_uses_RF .and. SIM_uses_RF) then
994 >       if ((FF_uses_RF .and. SIM_uses_RF) .or. (corrMethod == 3)) then
995  
996   #ifdef IS_MPI
997            call scatter(rf_Row,rf,plan_atom_row_3d)
# Line 1068 | Line 1075 | contains
1075      real ( kind = dp ), intent(inout) :: rijsq
1076      real ( kind = dp )                :: r
1077      real ( kind = dp ), intent(inout) :: d(3)
1071    real ( kind = dp ) :: ebalance
1078      integer :: me_i, me_j
1079  
1080      integer :: iHash
# Line 1093 | Line 1099 | contains
1099  
1100      if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1101         call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, &
1102 <            pot, eFrame, f, t, do_pot, corrMethod)
1102 >            pot, eFrame, f, t, do_pot, corrMethod, rcuti)
1103  
1104 <       if (FF_uses_RF .and. SIM_uses_RF) then
1104 >       if ((FF_uses_RF .and. SIM_uses_RF) .or. (corrMethod == 3)) then
1105  
1106            ! CHECK ME (RF needs to know about all electrostatic types)
1107            call accumulate_rf(i, j, r, eFrame, sw)
# Line 1378 | Line 1384 | contains
1384    function FF_RequiresPostpairCalc() result(doesit)
1385      logical :: doesit
1386      doesit = FF_uses_RF
1387 +    if (corrMethod == 3) doesit = .true.
1388    end function FF_RequiresPostpairCalc
1389  
1390   #ifdef PROFILE

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines