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

Comparing trunk/OOPSE-4/src/UseTheForce/doForces.F90 (file contents):
Revision 2280 by gezelter, Thu Sep 1 20:17:55 2005 UTC vs.
Revision 2281 by gezelter, Thu Sep 1 22:56:20 2005 UTC

# Line 45 | Line 45
45  
46   !! @author Charles F. Vardeman II
47   !! @author Matthew Meineke
48 < !! @version $Id: doForces.F90,v 1.34 2005-09-01 20:17:55 gezelter Exp $, $Date: 2005-09-01 20:17:55 $, $Name: not supported by cvs2svn $, $Revision: 1.34 $
48 > !! @version $Id: doForces.F90,v 1.35 2005-09-01 22:56:20 gezelter Exp $, $Date: 2005-09-01 22:56:20 $, $Name: not supported by cvs2svn $, $Revision: 1.35 $
49  
50  
51   module doForces
# Line 275 | Line 275 | contains
275      nAtypes = getSize(atypes)
276      
277      do i = 1, nAtypes
278 <       if (SimHasAtype(i)) then          
278 >       if (SimHasAtype(i)) then    
279            call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
280            call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
281            call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
# Line 284 | Line 284 | contains
284            call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
285            call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
286            
287 +          atypeMaxCutoff(i) = 0.0_dp
288            if (i_is_LJ) then
289               thisRcut = getSigma(i) * 2.5_dp
290               if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
# Line 327 | Line 328 | contains
328   #else
329      iend = nGroups
330   #endif
331 <
331 >    
332      !! allocate the groupToGtype and gtypeMaxCutoff here.
333      
334 <    !! first we do a single loop over the cutoff groups to find the largest cutoff for any atypes
335 <    !! present in this group.   We also create gtypes at this point.
334 >    !! first we do a single loop over the cutoff groups to find the
335 >    !! largest cutoff for any atypes present in this group.  We also
336 >    !! create gtypes at this point.
337 >    
338      tol = 1.0d-6
339      
340      do i = istart, iend      
# Line 366 | Line 369 | contains
369      enddo
370      
371      !! allocate the gtypeCutoffMap here.
372 <
372 >    
373      !! then we do a double loop over all the group TYPES to find the cutoff
374      !! map between groups of two types
375      
# Line 374 | Line 377 | contains
377         do j = 1, nGroupTypes
378        
379            select case(cutoffPolicy)
380 <             case(TRADITIONAL_CUTOFF_POLICY)
381 <                thisRcut = maxval(gtypeMaxCutoff)
382 <             case(MIX_CUTOFF_POLICY)
383 <                thisRcut = 0.5_dp * (gtypeMaxCutoff(i) + gtypeMaxCutoff(j))
384 <             case(MAX_CUTOFF_POLICY)
385 <                thisRcut = max(gtypeMaxCutoff(i), gtypeMaxCutoff(j))
386 <             case default
387 <                call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
388 <                return
389 <          end select      
390 <         gtypeCutoffMap(i,j)%rcut = thisRcut
391 <         gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
392 <         skin = defaultRlist - defaultRcut
393 <         gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skin)**2
380 >          case(TRADITIONAL_CUTOFF_POLICY)
381 >             thisRcut = maxval(gtypeMaxCutoff)
382 >          case(MIX_CUTOFF_POLICY)
383 >             thisRcut = 0.5_dp * (gtypeMaxCutoff(i) + gtypeMaxCutoff(j))
384 >          case(MAX_CUTOFF_POLICY)
385 >             thisRcut = max(gtypeMaxCutoff(i), gtypeMaxCutoff(j))
386 >          case default
387 >             call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
388 >             return
389 >          end select
390 >          gtypeCutoffMap(i,j)%rcut = thisRcut
391 >          gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
392 >          skin = defaultRlist - defaultRcut
393 >          gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skin)**2
394         enddo
395      enddo
396      
397      haveGtypeCutoffMap = .true.
398 <   end subroutine createGtypeCutoffMap
399 <
400 <   subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
401 <     real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
402 <     integer, intent(in) :: cutPolicy
403 <
404 <     defaultRcut = defRcut
405 <     defaultRsw = defRsw
406 <     defaultRlist = defRlist
407 <     cutoffPolicy = cutPolicy
408 <   end subroutine setDefaultCutoffs
409 <
410 <   subroutine setCutoffPolicy(cutPolicy)
398 >  end subroutine createGtypeCutoffMap
399 >  
400 >  subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
401 >    real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
402 >    integer, intent(in) :: cutPolicy
403 >    
404 >    defaultRcut = defRcut
405 >    defaultRsw = defRsw
406 >    defaultRlist = defRlist
407 >    cutoffPolicy = cutPolicy
408 >  end subroutine setDefaultCutoffs
409 >  
410 >  subroutine setCutoffPolicy(cutPolicy)
411  
412       integer, intent(in) :: cutPolicy
413       cutoffPolicy = cutPolicy

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines