FEI Version of the Day
fei_EqnComm.cpp
00001 /*
00002 // @HEADER
00003 // ************************************************************************
00004 //             FEI: Finite Element Interface to Linear Solvers
00005 //                  Copyright (2005) Sandia Corporation.
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, the
00008 // U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Alan Williams (william@sandia.gov) 
00038 //
00039 // ************************************************************************
00040 // @HEADER
00041 */
00042 
00043 
00044 #include "fei_EqnComm.hpp"
00045 #include "fei_sstream.hpp"
00046 
00047 namespace fei {
00048 
00049 EqnComm::EqnComm(MPI_Comm comm, int numLocalEqns)
00050  : comm_(comm),
00051    globalOffsets_(2)
00052 {
00053 #ifndef FEI_SER
00054 
00055   int numProcs, localProc;
00056   MPI_Comm_size(comm, &numProcs);
00057   MPI_Comm_rank(comm, &localProc);
00058 
00059   std::vector<int> local(numProcs*2, 0);
00060   int* global = &local[0] + numProcs;
00061 
00062   if (numLocalEqns < 0) {
00063     throw std::runtime_error("fei::EqnComm ERROR, negative numLocalEqns not allowed.");
00064   }
00065 
00066   local[localProc] = numLocalEqns;
00067 
00068   MPI_Allreduce(&local[0], global, numProcs, MPI_INT, MPI_MAX, comm_);
00069 
00070   globalOffsets_.resize(numProcs+1);
00071 
00072   int offset = 0;
00073   for(int i=0; i<numProcs; ++i) {
00074     globalOffsets_[i] = offset;
00075     offset += global[i];
00076   }
00077   globalOffsets_[numProcs] = offset;
00078 
00079 #else
00080   globalOffsets_[0] = 0;
00081   globalOffsets_[1] = numLocalEqns;
00082 #endif
00083 }
00084   
00085 EqnComm::EqnComm(MPI_Comm comm, int numLocalEqns, const std::vector<int>& globalOffsets)
00086  : comm_(comm),
00087    globalOffsets_(globalOffsets)
00088 {
00089 }
00090   
00091 EqnComm::~EqnComm()
00092 {
00093 }
00094 
00095 const std::vector<int>&
00096 EqnComm::getGlobalOffsets() const
00097 {
00098   return(globalOffsets_);
00099 }
00100 
00101 int
00102 EqnComm::getOwnerProc(int eqn) const
00103 {
00104 //  std::vector<int>::const_iterator
00105 //   iter = std::lower_bound(globalOffsets_.begin(), globalOffsets_.end(),
00106 //                           eqn);
00107 //  int proc = iter - globalOffsets_.begin() - 1;
00108 //  if (*iter==eqn) ++proc;
00109   int proc = 0;
00110   for(unsigned p=1; p<globalOffsets_.size(); ++p) {
00111     if (eqn < globalOffsets_[p]) {
00112       proc = p-1;
00113       break;
00114     }
00115   }
00116 
00117 #ifndef NDEBUG
00118   int numProcs = globalOffsets_.size()-1;
00119   if (proc >= numProcs) {
00120     FEI_OSTRINGSTREAM osstr;
00121     osstr << "fei::EqnComm::getOwnerProc: input eqn="<<eqn<<", proc="<<proc
00122       << ", ERROR, proc should be in [0.."<<numProcs-1<<"].";
00123     throw std::runtime_error(std::string(osstr.str().c_str()));
00124   }
00125 #endif
00126 
00127   return((int)proc);
00128 }
00129 
00130 }//namespace fei
00131 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Friends