ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/do_Forces.F90 (file contents):
Revision 894 by chuckv, Mon Jan 5 21:00:05 2004 UTC vs.
Revision 895 by chuckv, Mon Jan 5 22:12:11 2004 UTC

# Line 4 | Line 4
4  
5   !! @author Charles F. Vardeman II
6   !! @author Matthew Meineke
7 < !! @version $Id: do_Forces.F90,v 1.41 2004-01-05 21:00:05 chuckv Exp $, $Date: 2004-01-05 21:00:05 $, $Name: not supported by cvs2svn $, $Revision: 1.41 $
7 > !! @version $Id: do_Forces.F90,v 1.42 2004-01-05 22:12:11 chuckv Exp $, $Date: 2004-01-05 22:12:11 $, $Name: not supported by cvs2svn $, $Revision: 1.42 $
8  
9   module do_Forces
10    use force_globals
# Line 51 | Line 51 | module do_Forces
51    real :: forceTimeInitial, forceTimeFinal
52    integer :: nLoops
53   #endif
54 +
55 +  logical, allocatable :: propertyMapI(:,:)
56 +  logical, allocatable :: propertyMapJ(:,:)
57  
58   contains
59  
# Line 204 | Line 207 | contains
207      real ( kind = dp ), dimension(3,getNlocal()) :: f
208      !! Torsion array provided by C, dimensioned by getNlocal
209      real( kind = dp ), dimension(3,getNlocal()) :: t    
210 +
211      !! Stress Tensor
212      real( kind = dp), dimension(9) :: tau  
213      real ( kind = dp ) :: pot
# Line 255 | Line 259 | contains
259      do_pot = do_pot_c
260      do_stress = do_stress_c
261  
262 <
262 >
263 > #ifdef IS_MPI
264 >    if (.not.allocated(propertyMapI)) then
265 >       allocate(propertyMapI(5,getNrow())
266 >    endif
267 >
268 >    do i = 1, nrow
269 >       me_i = atid_row(i)
270 > #else
271 >    if (.not.allocated(propertyMapI)) then
272 >       allocate(propertyMapI(5,getNlocal())
273 >    endif
274 >
275 >    do i = 1, natoms
276 >       me_i = atid(i)
277 > #endif
278 >      
279 >       propertyMapI(1:5,i) = .false.
280 >
281 >       call getElementProperty(atypes, me_i, "propertyPack", propPack_i)
282 >    
283 >       ! unpack the properties
284 >      
285 >       if (iand(propPack_i, LJ_PROPERTY_MASK) .eq. LJ_PROPERTY_MASK) &
286 >            propertyMapI(1, i) = .true.
287 >       if (iand(propPack_i, DP_PROPERTY_MASK) .eq. DP_PROPERTY_MASK) &
288 >            propertyMapI(2, i) = .true.
289 >       if (iand(propPack_i, STICKY_PROPERTY_MASK) .eq. STICKY_PROPERTY_MASK) &
290 >            propertyMapI(3, i) = .true.
291 >       if (iand(propPack_i, GB_PROPERTY_MASK) .eq. GB_PROPERTY_MASK) &
292 >            propertyMapI(4, i) = .true.
293 >       if (iand(propPack_i, EAM_PROPERTY_MASK) .eq. EAM_PROPERTY_MASK) &
294 >            propertyMapI(5, i) = .true.
295 >
296 >    end do
297 >
298 > #ifdef IS_MPI
299 >    if (.not.allocated(propertyMapJ)) then
300 >       allocate(propertyMapJ(5,getNcol())
301 >    endif
302 >
303 >    do j = 1, ncol
304 >       me_j = atid_col(j)
305 > #else
306 >    if (.not.allocated(propertyMapJ)) then
307 >       allocate(propertyMapJ(5,getNlocal())
308 >    endif
309 >
310 >    do j = 1, natoms
311 >       me_j = atid(j)
312 > #endif
313 >      
314 >       propertyMapJ(1:5,j) = .false.
315 >
316 >       call getElementProperty(atypes, me_j, "propertyPack", propPack_j)
317 >    
318 >       ! unpack the properties
319 >      
320 >       if (iand(propPack_j, LJ_PROPERTY_MASK) .eq. LJ_PROPERTY_MASK) &
321 >            propertyMapJ(1, j) = .true.
322 >       if (iand(propPack_j, DP_PROPERTY_MASK) .eq. DP_PROPERTY_MASK) &
323 >            propertyMapJ(2, j) = .true.
324 >       if (iand(propPack_j, STICKY_PROPERTY_MASK) .eq. STICKY_PROPERTY_MASK) &
325 >            propertyMapJ(3, j) = .true.
326 >       if (iand(propPack_j, GB_PROPERTY_MASK) .eq. GB_PROPERTY_MASK) &
327 >            propertyMapJ(4, j) = .true.
328 >       if (iand(propPack_j, EAM_PROPERTY_MASK) .eq. EAM_PROPERTY_MASK) &
329 >            propertyMapJ(5, j) = .true.
330 >
331 >    end do
332 >
333      ! Gather all information needed by all force loops:
334      
335   #ifdef IS_MPI    
# Line 451 | Line 525 | contains
525         nlist = 0      
526        
527         do i = 1, nrow
528 +
529            point(i) = nlist + 1
530            
531            inner: do j = 1, ncol
# Line 739 | Line 814 | contains
814      me_j = atid(j)
815  
816   #endif
742
743    call getElementProperty(atypes, me_i, "propertyPack", propPack_i)
744    call getElementProperty(atypes, me_j, "propertyPack", propPack_j)
817      
746    ! unpack the properties
747
748    if (iand(propPack_i, LJ_PROPERTY_MASK) .eq. LJ_PROPERTY_MASK) is_LJ_i = .true.
749    if (iand(propPack_i, DP_PROPERTY_MASK) .eq. DP_PROPERTY_MASK) is_DP_i = .true.
750    if (iand(propPack_i, STICKY_PROPERTY_MASK) .eq. STICKY_PROPERTY_MASK) is_Sticky_i = .true.
751    if (iand(propPack_i, GB_PROPERTY_MASK) .eq. GB_PROPERTY_MASK) is_GB_i = .true.
752    if (iand(propPack_i, EAM_PROPERTY_MASK) .eq. EAM_PROPERTY_MASK) is_EAM_i = .true.
753
754    if (iand(propPack_j, LJ_PROPERTY_MASK) .eq. LJ_PROPERTY_MASK) is_LJ_j = .true.
755    if (iand(propPack_j, DP_PROPERTY_MASK) .eq. DP_PROPERTY_MASK) is_DP_j = .true.
756    if (iand(propPack_j, STICKY_PROPERTY_MASK) .eq. STICKY_PROPERTY_MASK) is_Sticky_j = .true.
757    if (iand(propPack_j, GB_PROPERTY_MASK) .eq. GB_PROPERTY_MASK) is_GB_j = .true.
758    if (iand(propPack_j, EAM_PROPERTY_MASK) .eq. EAM_PROPERTY_MASK) is_EAM_j = .true.
759
760
818      if (FF_uses_LJ .and. SimUsesLJ()) then
819  
820 <       if ( is_LJ_i .and. is_LJ_j ) &
820 >       if ( propertyMapI(1, me_i) .and. propertyMapJ(1, me_j) ) &
821              call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
822  
823      endif
824  
825      if (FF_uses_dipoles .and. SimUsesDipoles()) then
826        
827 <       if ( is_DP_i .and. is_DP_j ) then
827 >       if ( propertyMapI(2, me_i) .and. propertyMapJ(2, me_j)) then
828            call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
829                 do_pot, do_stress)
830            if (FF_uses_RF .and. SimUsesRF()) then
# Line 780 | Line 837 | contains
837  
838      if (FF_uses_Sticky .and. SimUsesSticky()) then
839  
840 <       if ( is_Sticky_i .and. is_Sticky_j ) then
840 >       if ( propertyMapI(3, me_i) .and. propertyMapJ(3, me_j)) then
841            call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
842                 do_pot, do_stress)
843         endif
# Line 789 | Line 846 | contains
846  
847      if (FF_uses_GB .and. SimUsesGB()) then
848        
849 <       if ( is_GB_i .and. is_GB_j ) then
849 >       if ( propertyMapI(4, me_i) .and. propertyMapJ(4, me_j)) then
850            call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
851                 do_pot, do_stress)          
852         endif
# Line 800 | Line 857 | contains
857    
858     if (FF_uses_EAM .and. SimUsesEAM()) then
859        
860 <      if ( is_EAM_i .and. is_EAM_j ) &
861 <           call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
860 >      if ( propertyMapI(5, me_i) .and. propertyMapJ(5, me_j)) then
861 >         call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
862 >      endif
863  
864     endif
865  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines