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

Comparing trunk/OOPSE/libmdtools/mpiSimulation_module.F90 (file contents):
Revision 1203 by gezelter, Thu May 27 18:59:17 2004 UTC vs.
Revision 1214 by gezelter, Tue Jun 1 18:42:58 2004 UTC

# Line 7 | Line 7
7   !!
8   !! @author Charles F. Vardeman II
9   !! @author Matthew Meineke
10 < !! @version $Id: mpiSimulation_module.F90,v 1.14 2004-05-27 18:59:15 gezelter Exp $, $Date: 2004-05-27 18:59:15 $, $Name: not supported by cvs2svn $, $Revision: 1.14 $
10 > !! @version $Id: mpiSimulation_module.F90,v 1.15 2004-06-01 18:42:58 gezelter Exp $, $Date: 2004-06-01 18:42:58 $, $Name: not supported by cvs2svn $, $Revision: 1.15 $
11  
12   module mpiSimulation  
13    use definitions
# Line 141 | Line 141 | contains
141  
142    !! Sets up mpiComponentPlan with structure passed from C++.
143    subroutine setupSimParallel(thisComponentPlan, nAtomTags, atomTags, &
144 <       status)
144 >       nGroupTags, groupTags, status)
145      !! Passed Arguments
146      !! mpiComponentPlan struct from C
147      type (mpiComponentPlan), intent(inout) :: thisComponentPlan
148      !! Number of tags passed
149 <    integer, intent(in) :: nAtomTags
149 >    integer, intent(in) :: nAtomTags, nGroupTags
150      !! Result status, 0 = normal, -1 = error
151      integer, intent(out) :: status
152      integer :: localStatus
153      !! Global reference tag for local particles
154      integer, dimension(nAtomTags), intent(inout) :: atomTags
155 +    integer, dimension(nGroupTags), intent(inout) :: groupTags
156  
157      !write(*,*) 'mpiSim_mod thinks node', thisComponentPlan%myNode, &
158      !     ' has atomTags(1) = ', atomTags(1)
# Line 210 | Line 211 | contains
211         status = -1
212         return
213      endif
214 +
215 +
216 +    call setGroupTags(groupTags,localStatus)
217 +    if (localStatus /= 0) then
218 +       status = -1
219 +       return
220 +    endif
221 +
222      isSimSet = .true.
223  
224   !    call printComponentPlan(mpiSim,0)
# Line 708 | Line 717 | contains
717      call gather(tags, AtomColToGlobal, plan_atom_col)
718      
719    end subroutine setAtomTags
720 +
721 +  subroutine setGroupTags(tags, status)
722 +    integer, dimension(:) :: tags
723 +    integer :: status
724 +
725 +    integer :: alloc_stat
726 +    
727 +    integer :: nGroupsInCol
728 +    integer :: nGroupsInRow
729 +
730 +    status = 0
731 +
732 +    nGroupsInRow = getNgroupsInRow(plan_group_row)
733 +    nGroupsInCol = getNgroupsInCol(plan_group_col)
734 +    
735 +    if(allocated(GroupLocalToGlobal)) then
736 +       deallocate(GroupLocalToGlobal)
737 +    endif
738 +    allocate(GroupLocalToGlobal(size(tags)),STAT=alloc_stat)
739 +    if (alloc_stat /= 0 ) then
740 +       status = -1
741 +       return
742 +    endif
743 +    
744 +    GroupLocalToGlobal = tags
745 +
746 +    if(allocated(GroupRowToGlobal)) then
747 +       deallocate(GroupRowToGlobal)
748 +    endif
749 +    allocate(GroupRowToGlobal(nGroupsInRow),STAT=alloc_stat)
750 +    if (alloc_stat /= 0 ) then
751 +       status = -1
752 +       return
753 +    endif
754 +
755 +    if(allocated(GroupColToGlobal)) then
756 +       deallocate(GroupColToGlobal)
757 +    endif
758 +    allocate(GroupColToGlobal(nGroupsInCol),STAT=alloc_stat)
759 +    if (alloc_stat /= 0 ) then
760 +       status = -1
761 +       return
762 +    endif
763 +    
764 +    call gather(tags, GroupRowToGlobal, plan_group_row)
765 +    call gather(tags, GroupColToGlobal, plan_group_col)
766 +    
767 +  end subroutine setGroupTags
768    
769    function getNatomsInCol(thisplan) result(nInCol)
770      type (gs_plan), intent(in) :: thisplan

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines