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 2279 by chrisfen, Tue Aug 30 18:23:50 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.33 2005-08-30 18:23:29 chrisfen Exp $, $Date: 2005-08-30 18:23:29 $, $Name: not supported by cvs2svn $, $Revision: 1.33 $
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 85 | Line 85 | module doForces
85    logical, save :: haveSaneForceField = .false.
86    logical, save :: haveInteractionHash = .false.
87    logical, save :: haveGtypeCutoffMap = .false.
88 +  logical, save :: haveRlist = .false.
89  
90    logical, save :: FF_uses_DirectionalAtoms
91    logical, save :: FF_uses_Dipoles
# Line 257 | Line 258 | contains
258      logical :: i_is_Shape
259  
260      integer :: myStatus, nAtypes,  i, j, istart, iend, jstart, jend
261 <    integer :: n_in_i
262 <    real(kind=dp):: thisSigma, bigSigma, thisRcut
261 >    integer :: n_in_i, me_i, ia, g, atom1, nGroupTypes
262 >    real(kind=dp):: thisSigma, bigSigma, thisRcut, tol, skin
263      real(kind=dp) :: biggestAtypeCutoff
264  
265      stat = 0
# Line 274 | 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 283 | 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 317 | Line 319 | contains
319            endif
320         endif
321      enddo
322 <
322 >  
323 >    nGroupTypes = 0
324 >    
325      istart = 1
326   #ifdef IS_MPI
327      iend = nGroupsInRow
328   #else
329      iend = nGroups
330   #endif
331 <    outer: do i = istart, iend
332 <      
333 <       n_in_i = groupStartRow(i+1) - groupStartRow(i)
334 <      
331 >    
332 >    !! allocate the groupToGtype and gtypeMaxCutoff here.
333 >    
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      
341 >       n_in_i = groupStartRow(i+1) - groupStartRow(i)
342 >       groupMaxCutoff(i) = 0.0_dp
343 >       do ia = groupStartRow(i), groupStartRow(i+1)-1
344 >          atom1 = groupListRow(ia)
345   #ifdef IS_MPI
346 <       jstart = 1
333 <       jend = nGroupsInCol
346 >          me_i = atid_row(atom1)
347   #else
348 <       jstart = i+1
349 <       jend = nGroups
350 < #endif
351 <      
352 <      
353 <      
354 <      
355 <      
356 <      
357 <    enddo outer        
348 >          me_i = atid(atom1)
349 > #endif          
350 >          if (atypeMaxCutoff(me_i).gt.groupMaxCutoff(i)) then
351 >             groupMaxCutoff(i)=atypeMaxCutoff(me_i)
352 >          endif
353 >       enddo
354 >       if (nGroupTypes.eq.0) then
355 >          nGroupTypes = nGroupTypes + 1
356 >          gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
357 >          groupToGtype(i) = nGroupTypes
358 >       else
359 >          do g = 1, nGroupTypes
360 >             if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).gt.tol) then
361 >                nGroupTypes = nGroupTypes + 1
362 >                gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
363 >                groupToGtype(i) = nGroupTypes
364 >             else
365 >                groupToGtype(i) = g
366 >             endif
367 >          enddo
368 >       endif
369 >    enddo
370      
371 <     haveGtypeCutoffMap = .true.
372 <   end subroutine createGtypeCutoffMap
373 <
374 <   subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy)
375 <     real(kind=dp),intent(in) :: defRcut, defRsw, defRlist
376 <     integer, intent(in) :: cutPolicy
377 <
378 <     defaultRcut = defRcut
379 <     defaultRsw = defRsw
380 <     defaultRlist = defRlist
381 <     cutoffPolicy = cutPolicy
382 <   end subroutine setDefaultCutoffs
383 <
384 <   subroutine setCutoffPolicy(cutPolicy)
371 >    !! allocate the gtypeCutoffMap here.
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 >    
376 >    do i = 1, nGroupTypes
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
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)
411  
412       integer, intent(in) :: cutPolicy
413       cutoffPolicy = cutPolicy
# Line 438 | Line 489 | contains
489    end subroutine doReadyCheck
490  
491  
492 <  subroutine init_FF(use_RF_c, use_UW_c, use_DW_c, thisStat)
492 >  subroutine init_FF(use_RF, use_UW, use_DW, thisStat)
493  
494 <    logical, intent(in) :: use_RF_c
495 <    logical, intent(in) :: use_UW_c
496 <    logical, intent(in) :: use_DW_c
494 >    logical, intent(in) :: use_RF
495 >    logical, intent(in) :: use_UW
496 >    logical, intent(in) :: use_DW
497      integer, intent(out) :: thisStat  
498      integer :: my_status, nMatches
499      integer :: corrMethod
# Line 453 | Line 504 | contains
504      thisStat = 0
505  
506      !! Fortran's version of a cast:
507 <    FF_uses_RF = use_RF_c
507 >    FF_uses_RF = use_RF
508  
509      !! set the electrostatic correction method
510 <    if (use_UW_c .eq. .true.) then
510 >    if (use_UW) then
511         corrMethod = 1
512 <    elseif (use_DW_c .eq. .true.) then
512 >    elseif (use_DW) then
513         corrMethod = 2
514      else
515         corrMethod = 0
# Line 719 | Line 770 | contains
770                    q_group(:,j), d_grp, rgrpsq)
771   #endif
772  
773 <             if (rgrpsq < InteractionHash(me_i,me_j)%rListsq) then
773 >             if (rgrpsq < gtypeCutoffMap(groupToGtype(i),groupToGtype(j))%rListsq) then
774                  if (update_nlist) then
775                     nlist = nlist + 1
776  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines