EpetraExt_MultiMpiComm.cpp

Go to the documentation of this file.
00001 //@HEADER
00002 /*
00003 ************************************************************************
00004 
00005               EpetraExt: Extended Linear Algebra Services Package 
00006                 Copyright (2001) Sandia Corporation
00007 
00008 Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00009 license for use of this work by or on behalf of the U.S. Government.
00010 
00011 This library is free software; you can redistribute it and/or modify
00012 it under the terms of the GNU Lesser General Public License as
00013 published by the Free Software Foundation; either version 2.1 of the
00014 License, or (at your option) any later version.
00015  
00016 This library is distributed in the hope that it will be useful, but
00017 WITHOUT ANY WARRANTY; without even the implied warranty of
00018 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00019 Lesser General Public License for more details.
00020  
00021 You should have received a copy of the GNU Lesser General Public
00022 License along with this library; if not, write to the Free Software
00023 Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00024 USA
00025 Questions? Contact Michael A. Heroux (maherou@sandia.gov) 
00026 
00027 ************************************************************************
00028 */
00029 //@HEADER
00030 
00031 #include "EpetraExt_MultiMpiComm.h" 
00032 
00033 namespace EpetraExt {
00034 
00035 MultiMpiComm::MultiMpiComm(MPI_Comm globalMpiComm, int subDomainProcs, int numTimeSteps_) :
00036         Epetra_MpiComm(globalMpiComm), subComm(0), numSubDomains(-1),
00037         subDomainRank(-1), numTimeSteps(-1),
00038   numTimeStepsOnDomain(-1), firstTimeStepOnDomain(-1)
00039 {
00040   //Need to construct subComm for each sub domain, compute subDomainRank,
00041   //and check that all integer arithmatic works out correctly.
00042  
00043   int ierrmpi, size, rank;
00044   ierrmpi = MPI_Comm_size(globalMpiComm, &size);
00045   ierrmpi = MPI_Comm_rank(globalMpiComm, &rank);
00046 
00047   if (size % subDomainProcs != 0) {cout<<"ERROR: num subDomainProcs "<< subDomainProcs
00048      << " does not divide into num total procs " << size << endl; exit(-1);}
00049 
00050   numSubDomains = size / subDomainProcs;
00051 
00052   // Create split communicators, the size of subDomainProcs
00053   MPI_Comm split_MPI_Comm;
00054   subDomainRank = rank/subDomainProcs;
00055   ierrmpi =  MPI_Comm_split(globalMpiComm, subDomainRank, rank, &split_MPI_Comm);
00056 
00057   // Construct second epetra communicators
00058   subComm = new Epetra_MpiComm(split_MPI_Comm);
00059 
00060   // Compute number of time steps on this sub domain
00061   ResetNumTimeSteps(numTimeSteps_);
00062 
00063   if (numTimeSteps_ > 0)
00064     cout << "Processor " << rank << " is on subdomain " << subDomainRank 
00065          << " and owns " << numTimeStepsOnDomain << " time steps, starting with " 
00066          <<  firstTimeStepOnDomain << endl;
00067   else
00068     cout << "Processor " << rank << " is on subdomain " << subDomainRank << endl;
00069 }
00070 
00071 // This constructor is for just one subdomain, so only adds the info
00072 // for multiple time steps on the domain. No two-level parallelism.
00073 MultiMpiComm::MultiMpiComm(const Epetra_MpiComm& EpetraMpiComm_, int numTimeSteps_) :
00074         Epetra_MpiComm(EpetraMpiComm_), subComm(0), numSubDomains(1),
00075         subDomainRank(0), numTimeSteps(numTimeSteps_),
00076   numTimeStepsOnDomain(numTimeSteps_), firstTimeStepOnDomain(0)
00077 {
00078    subComm = new Epetra_MpiComm(EpetraMpiComm_);
00079 }
00080   
00081 //Copy Constructor
00082 MultiMpiComm::MultiMpiComm(const MultiMpiComm &MMC ) :
00083         Epetra_MpiComm(MMC), subComm(new Epetra_MpiComm(*(MMC.subComm))),
00084         numSubDomains(MMC.numSubDomains), subDomainRank(MMC.subDomainRank),
00085   numTimeSteps(MMC.numTimeSteps), 
00086   numTimeStepsOnDomain(MMC.numTimeStepsOnDomain), 
00087   firstTimeStepOnDomain(MMC.firstTimeStepOnDomain) 
00088 {
00089 }
00090 
00091 MultiMpiComm::~MultiMpiComm()
00092 {
00093   delete subComm;
00094 }
00095 
00096 void MultiMpiComm::ResetNumTimeSteps(int numTimeSteps_)
00097 {
00098   numTimeSteps = numTimeSteps_;
00099 
00100   // Compute number of time steps on this sub domain
00101   if (numTimeSteps > 0) {
00102     // Compute part for number of domains dividing evenly into number of steps
00103     numTimeStepsOnDomain = numTimeSteps / numSubDomains; 
00104     firstTimeStepOnDomain = numTimeStepsOnDomain * subDomainRank;
00105 
00106     // Dole out remainder
00107     int remainder = numTimeSteps % numSubDomains;
00108     if (subDomainRank < remainder) {
00109       numTimeStepsOnDomain++; 
00110       firstTimeStepOnDomain += subDomainRank; 
00111     }
00112     else firstTimeStepOnDomain += remainder; 
00113   }
00114   else {
00115     numTimeStepsOnDomain = -1;
00116     firstTimeStepOnDomain = -1;
00117   }
00118 }
00119 
00120 } //namespace EpetraExt

Generated on Tue Oct 20 12:45:30 2009 for EpetraExt by doxygen 1.4.7