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 2342 by chrisfen, Tue Oct 4 19:32:58 2005 UTC vs.
Revision 2350 by chuckv, Mon Oct 10 21:20:46 2005 UTC

# Line 45 | Line 45
45  
46   !! @author Charles F. Vardeman II
47   !! @author Matthew Meineke
48 < !! @version $Id: doForces.F90,v 1.50 2005-10-04 19:32:58 chrisfen Exp $, $Date: 2005-10-04 19:32:58 $, $Name: not supported by cvs2svn $, $Revision: 1.50 $
48 > !! @version $Id: doForces.F90,v 1.51 2005-10-10 21:20:46 chuckv Exp $, $Date: 2005-10-10 21:20:46 $, $Name: not supported by cvs2svn $, $Revision: 1.51 $
49  
50  
51   module doForces
# Line 124 | Line 124 | module doForces
124    ! Bit hash to determine pair-pair interactions.
125    integer, dimension(:,:), allocatable :: InteractionHash
126    real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff
127 <  real(kind=dp), dimension(:), allocatable :: groupMaxCutoff
128 <  integer, dimension(:), allocatable :: groupToGtype
129 <  real(kind=dp), dimension(:), allocatable :: gtypeMaxCutoff
127 >  real(kind=dp), dimension(:), allocatable, target :: groupMaxCutoffRow
128 >  real(kind=dp), dimension(:), pointer :: groupMaxCutoffCol
129 >
130 >  integer, dimension(:), allocatable, target :: groupToGtypeRow
131 >  integer, dimension(:), pointer :: groupToGtypeCol => null()
132 >
133 >  real(kind=dp), dimension(:), allocatable,target :: gtypeMaxCutoffRow
134 >  real(kind=dp), dimension(:), pointer :: gtypeMaxCutoffCol
135    type ::gtypeCutoffs
136       real(kind=dp) :: rcut
137       real(kind=dp) :: rcutsq
# Line 260 | Line 265 | contains
265      logical :: GtypeFound
266  
267      integer :: myStatus, nAtypes,  i, j, istart, iend, jstart, jend
268 <    integer :: n_in_i, me_i, ia, g, atom1, nGroupTypes
268 >    integer :: n_in_i, me_i, ia, g, atom1
269      integer :: nGroupsInRow
270 <    real(kind=dp):: thisSigma, bigSigma, thisRcut, tol, skin
270 >    integer :: nGroupsInCol
271 >    integer :: nGroupTypesRow,nGroupTypesCol
272 >    real(kind=dp):: thisSigma, bigSigma, thisRcut, tradRcut, tol, skin
273      real(kind=dp) :: biggestAtypeCutoff
274  
275      stat = 0
# Line 276 | Line 283 | contains
283      endif
284   #ifdef IS_MPI
285      nGroupsInRow = getNgroupsInRow(plan_group_row)
286 +    nGroupsInCol = getNgroupsInCol(plan_group_col)
287   #endif
288      nAtypes = getSize(atypes)
289   ! Set all of the initial cutoffs to zero.
# Line 332 | Line 340 | contains
340         endif
341      enddo
342    
343 <    nGroupTypes = 0
343 >
344      
345      istart = 1
346 +    jstart = 1
347   #ifdef IS_MPI
348      iend = nGroupsInRow
349 +    jend = nGroupsInCol
350   #else
351      iend = nGroups
352 +    jend = nGroups
353   #endif
354      
355      !! allocate the groupToGtype and gtypeMaxCutoff here.
356 <    if(.not.allocated(groupToGtype)) then
357 <       allocate(groupToGtype(iend))
358 <       allocate(groupMaxCutoff(iend))
359 <       allocate(gtypeMaxCutoff(iend))
360 <       groupMaxCutoff = 0.0_dp
361 <       gtypeMaxCutoff = 0.0_dp
356 >    if(.not.allocated(groupToGtypeRow)) then
357 >     !  allocate(groupToGtype(iend))
358 >       allocate(groupToGtypeRow(iend))
359 >    else
360 >       deallocate(groupToGtypeRow)
361 >       allocate(groupToGtypeRow(iend))
362 >    endif
363 >    if(.not.allocated(groupMaxCutoffRow)) then
364 >       allocate(groupMaxCutoffRow(iend))
365 >    else
366 >       deallocate(groupMaxCutoffRow)
367 >       allocate(groupMaxCutoffRow(iend))
368 >    end if
369 >
370 >    if(.not.allocated(gtypeMaxCutoffRow)) then
371 >       allocate(gtypeMaxCutoffRow(iend))
372 >    else
373 >       deallocate(gtypeMaxCutoffRow)
374 >       allocate(gtypeMaxCutoffRow(iend))
375      endif
376 +
377 +
378 + #ifdef IS_MPI
379 +       ! We only allocate new storage if we are in MPI because Ncol /= Nrow
380 +    if(.not.allocated(groupToGtypeCol)) then
381 +       allocate(groupToGtypeCol(jend))
382 +    else
383 +       deallocate(groupToGtypeCol)
384 +       allocate(groupToGtypeCol(jend))
385 +    end if
386 +
387 +    if(.not.allocated(groupToGtypeCol)) then
388 +       allocate(groupToGtypeCol(jend))
389 +    else
390 +       deallocate(groupToGtypeCol)
391 +       allocate(groupToGtypeCol(jend))
392 +    end if
393 +    if(.not.allocated(gtypeMaxCutoffCol)) then
394 +       allocate(gtypeMaxCutoffCol(jend))
395 +    else
396 +       deallocate(gtypeMaxCutoffCol)      
397 +       allocate(gtypeMaxCutoffCol(jend))
398 +    end if
399 +
400 +       groupMaxCutoffCol = 0.0_dp
401 +       gtypeMaxCutoffCol = 0.0_dp
402 +
403 + #endif
404 +       groupMaxCutoffRow = 0.0_dp
405 +       gtypeMaxCutoffRow = 0.0_dp
406 +
407 +
408      !! first we do a single loop over the cutoff groups to find the
409      !! largest cutoff for any atypes present in this group.  We also
410      !! create gtypes at this point.
411      
412      tol = 1.0d-6
413 <    
413 >    nGroupTypesRow = 0
414 >
415      do i = istart, iend      
416         n_in_i = groupStartRow(i+1) - groupStartRow(i)
417 <       groupMaxCutoff(i) = 0.0_dp
417 >       groupMaxCutoffRow(i) = 0.0_dp
418         do ia = groupStartRow(i), groupStartRow(i+1)-1
419            atom1 = groupListRow(ia)
420   #ifdef IS_MPI
# Line 365 | Line 422 | contains
422   #else
423            me_i = atid(atom1)
424   #endif          
425 <          if (atypeMaxCutoff(me_i).gt.groupMaxCutoff(i)) then
426 <             groupMaxCutoff(i)=atypeMaxCutoff(me_i)
425 >          if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then
426 >             groupMaxCutoffRow(i)=atypeMaxCutoff(me_i)
427            endif          
428         enddo
429  
430 <       if (nGroupTypes.eq.0) then
431 <          nGroupTypes = nGroupTypes + 1
432 <          gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
433 <          groupToGtype(i) = nGroupTypes
430 >       if (nGroupTypesRow.eq.0) then
431 >          nGroupTypesRow = nGroupTypesRow + 1
432 >          gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
433 >          groupToGtypeRow(i) = nGroupTypesRow
434         else
435            GtypeFound = .false.
436 <          do g = 1, nGroupTypes
437 <             if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).lt.tol) then
438 <                groupToGtype(i) = g
436 >          do g = 1, nGroupTypesRow
437 >             if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then
438 >                groupToGtypeRow(i) = g
439                  GtypeFound = .true.
440               endif
441            enddo
442            if (.not.GtypeFound) then            
443 <             nGroupTypes = nGroupTypes + 1
444 <             gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i)
445 <             groupToGtype(i) = nGroupTypes
443 >             nGroupTypesRow = nGroupTypesRow + 1
444 >             gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
445 >             groupToGtypeRow(i) = nGroupTypesRow
446 >          endif
447 >       endif
448 >    enddo    
449 >
450 > #ifdef IS_MPI
451 >    do j = jstart, jend      
452 >       n_in_j = groupStartCol(j+1) - groupStartCol(j)
453 >       groupMaxCutoffCol(j) = 0.0_dp
454 >       do ja = groupStartCol(j), groupStartCol(j+1)-1
455 >          atom1 = groupListCol(ja)
456 >
457 >          me_j = atid_col(atom1)
458 >
459 >          if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then
460 >             groupMaxCutoffCol(j)=atypeMaxCutoff(me_j)
461 >          endif          
462 >       enddo
463 >
464 >       if (nGroupTypesCol.eq.0) then
465 >          nGroupTypesCol = nGroupTypesCol + 1
466 >          gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
467 >          groupToGtypeCol(j) = nGroupTypesCol
468 >       else
469 >          GtypeFound = .false.
470 >          do g = 1, nGroupTypesCol
471 >             if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then
472 >                groupToGtypeCol(j) = g
473 >                GtypeFound = .true.
474 >             endif
475 >          enddo
476 >          if (.not.GtypeFound) then            
477 >             nGroupTypesCol = nGroupTypesCol + 1
478 >             gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
479 >             groupToGtypeCol(j) = nGroupTypesCol
480            endif
481         endif
482      enddo    
483  
484 + #else
485 + ! Set pointers to information we just found
486 +    nGroupTypesCol = nGroupTypesRow
487 +    groupToGtypeCol => groupToGtypeRow
488 +    gtypeMaxCutoffCol => gtypeMaxCutoffRow
489 +    groupMaxCutoffCol => groupMaxCutoffRow
490 + #endif
491 +
492 +
493 +
494 +
495 +
496      !! allocate the gtypeCutoffMap here.
497 <    allocate(gtypeCutoffMap(nGroupTypes,nGroupTypes))
497 >    allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol))
498      !! then we do a double loop over all the group TYPES to find the cutoff
499      !! map between groups of two types
500 <    
501 <    do i = 1, nGroupTypes
502 <       do j = 1, nGroupTypes
500 >    tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol))
501 >
502 >    do i = 1, nGroupTypesRow
503 >       do j = 1, nGroupTypesCol
504        
505            select case(cutoffPolicy)
506            case(TRADITIONAL_CUTOFF_POLICY)
507 <             thisRcut = maxval(gtypeMaxCutoff)
507 >             thisRcut = tradRcut
508            case(MIX_CUTOFF_POLICY)
509 <             thisRcut = 0.5_dp * (gtypeMaxCutoff(i) + gtypeMaxCutoff(j))
509 >             thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j))
510            case(MAX_CUTOFF_POLICY)
511 <             thisRcut = max(gtypeMaxCutoff(i), gtypeMaxCutoff(j))
511 >             thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j))
512            case default
513               call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
514               return
# Line 424 | Line 528 | contains
528            endif
529         enddo
530      enddo
531 <
531 >    if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow)
532 >    if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow)
533 >    if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff)
534 > #ifdef IS_MPI
535 >    if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol)
536 >    if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol)
537 > #endif
538 >    groupMaxCutoffCol => null()
539 >    gtypeMaxCutoffCol => null()
540 >    
541      haveGtypeCutoffMap = .true.
542     end subroutine createGtypeCutoffMap
543  
# Line 788 | Line 901 | contains
901                    q_group(:,j), d_grp, rgrpsq)
902   #endif      
903  
904 <             if (rgrpsq < gtypeCutoffMap(groupToGtype(i),groupToGtype(j))%rListsq) then
904 >             if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
905                  if (update_nlist) then
906                     nlist = nlist + 1
907  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines