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

Comparing trunk/OOPSE-2.0/src/UseTheForce/doForces.F90 (file contents):
Revision 2285 by gezelter, Wed Sep 7 20:46:46 2005 UTC vs.
Revision 2290 by gezelter, Wed Sep 14 20:28:05 2005 UTC

# Line 45 | Line 45
45  
46   !! @author Charles F. Vardeman II
47   !! @author Matthew Meineke
48 < !! @version $Id: doForces.F90,v 1.38 2005-09-07 20:46:39 gezelter Exp $, $Date: 2005-09-07 20:46:39 $, $Name: not supported by cvs2svn $, $Revision: 1.38 $
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 $
49  
50  
51   module doForces
# Line 256 | Line 256 | contains
256      logical :: i_is_GB
257      logical :: i_is_EAM
258      logical :: i_is_Shape
259 +    logical :: GtypeFound
260  
261      integer :: myStatus, nAtypes,  i, j, istart, iend, jstart, jend
262      integer :: n_in_i, me_i, ia, g, atom1, nGroupTypes
263 +    integer :: nGroupsInRow
264      real(kind=dp):: thisSigma, bigSigma, thisRcut, tol, skin
265      real(kind=dp) :: biggestAtypeCutoff
266  
# Line 271 | Line 273 | contains
273            return
274         endif
275      endif
276 <
276 > #ifdef IS_MPI
277 >    nGroupsInRow = getNgroupsInRow(plan_group_row)
278 > #endif
279      nAtypes = getSize(atypes)
280      
281      do i = 1, nAtypes
# Line 353 | Line 357 | contains
357   #endif          
358            if (atypeMaxCutoff(me_i).gt.groupMaxCutoff(i)) then
359               groupMaxCutoff(i)=atypeMaxCutoff(me_i)
360 <          endif
360 >          endif          
361         enddo
362 +
363         if (nGroupTypes.eq.0) then
364            nGroupTypes = nGroupTypes + 1
365            gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
366            groupToGtype(i) = nGroupTypes
367         else
368 <          do g = 1, nGroupTypes
369 <             if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).gt.tol) then
370 <                nGroupTypes = nGroupTypes + 1
366 <                gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
367 <                groupToGtype(i) = nGroupTypes
368 <             else
368 >          GtypeFound = .false.
369 >          do g = 1, nGroupTypes
370 >             if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).lt.tol) then
371                  groupToGtype(i) = g
372 +                GtypeFound = .true.
373               endif
374            enddo
375 +          if (.not.GtypeFound) then            
376 +             nGroupTypes = nGroupTypes + 1
377 +             gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
378 +             groupToGtype(i) = nGroupTypes
379 +          endif
380         endif
381 <    enddo
382 <    
381 >    enddo    
382 >
383      !! allocate the gtypeCutoffMap here.
384      allocate(gtypeCutoffMap(nGroupTypes,nGroupTypes))
385      !! then we do a double loop over all the group TYPES to find the cutoff
# Line 380 | Line 388 | contains
388      do i = 1, nGroupTypes
389         do j = 1, nGroupTypes
390        
383          write(*,*) 'cutoffPolicy = ', cutoffPolicy
391            select case(cutoffPolicy)
392            case(TRADITIONAL_CUTOFF_POLICY)
393               thisRcut = maxval(gtypeMaxCutoff)
# Line 1151 | Line 1158 | contains
1158  
1159      integer :: me_i, me_j, iHash
1160  
1161 +    r = sqrt(rijsq)
1162 +
1163   #ifdef IS_MPI  
1164      me_i = atid_row(i)
1165      me_j = atid_col(j)  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines