FEI Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
EqnCommMgr Class Reference

#include <fei_EqnCommMgr.hpp>

Public Member Functions

 EqnCommMgr (MPI_Comm comm, bool accumulate=true)
 
 EqnCommMgr (const EqnCommMgr &src)
 
EqnCommMgroperator= (const EqnCommMgr &src)
 
virtual ~EqnCommMgr ()
 
EqnCommMgrdeepCopy ()
 
size_t getNumSharingProcs ()
 
std::vector< int > & sharingProcsPtr ()
 
size_t getNumOwnerProcs ()
 
std::vector< int > & ownerProcsPtr ()
 
void addLocalEqn (int eqnNumber, int srcProc)
 
void addSolnValues (int *eqnNumbers, double *values, int num)
 
int exchangeIndices (std::ostream *dbgOut=NULL)
 
int exchangeEqns (std::ostream *dbgOut=NULL)
 
void exchangeSoln ()
 
int mirrorProcEqns (ProcEqns &inProcEqns, ProcEqns &outProcEqns)
 
int mirrorProcEqnLengths (ProcEqns &inProcEqns, ProcEqns &outProcEqns)
 
int getNumLocalEqns ()
 
std::vector< int > & localEqnNumbers ()
 
std::vector< fei::CSVec * > & localEqns ()
 
std::vector< std::vector
< double > * > * 
localRHSsPtr ()
 
int addRemoteEqn (int eqnNumber, int destProc, const double *coefs, const int *indices, int num)
 
int addRemoteEqn (int eqnNumber, const double *coefs, const int *indices, int num)
 
int addRemoteEqns (fei::CSRMat &mat, bool onlyIndices)
 
int addRemoteRHS (fei::CSVec &vec, int rhsIndex)
 
void setNumRHSs (int numRHSs)
 
int addRemoteRHS (int eqnNumber, int destProc, int rhsIndex, double value)
 
int addRemoteRHS (int eqnNumber, int rhsIndex, double value)
 
void addRemoteIndices (int eqnNumber, int destProc, int *indices, int num)
 
int getNumRemoteEqns ()
 
std::vector< int > & sendEqnNumbersPtr ()
 
double * sendEqnSolnPtr ()
 
void resetCoefs ()
 
int gatherSharedBCs (EqnBuffer &bcEqns)
 
int exchangeRemEssBCs (int *essEqns, int numEssEqns, double *essAlpha, double *essGamma, MPI_Comm comm, std::ostream *dbgOut=NULL)
 
int getNumRemEssBCEqns ()
 
std::vector< int > & remEssBCEqnNumbersPtr ()
 
std::vector< fei::CSVec * > & remEssBCEqns ()
 
int exchangePtToBlkInfo (snl_fei::PointBlockMap &blkEqnMapper)
 
bool newCoefData ()
 
bool newRHSData ()
 
EqnBuffergetRecvEqns ()
 
EqnBuffergetSendEqns ()
 
ProcEqnsgetRecvProcEqns ()
 
ProcEqnsgetSendProcEqns ()
 

Static Public Member Functions

static int exchangeEqnBuffers (MPI_Comm comm, ProcEqns *sendProcEqns, EqnBuffer *sendEqns, ProcEqns *recvProcEqns, EqnBuffer *recvEqns, bool accumulate)
 

Public Attributes

bool accumulate_
 

Private Member Functions

void deleteEssBCs ()
 
int getSendProcNumber (int eqn)
 
int consistencyCheck (const char *caller, std::vector< int > &recvProcs, std::vector< int > &recvProcTotalLengths, std::vector< int > &sendProcs, std::vector< int > &sendProcTotalLengths)
 

Private Attributes

int localProc_
 
ProcEqnsrecvProcEqns_
 
bool exchangeIndicesCalled_
 
EqnBufferrecvEqns_
 
std::vector< double > solnValues_
 
ProcEqnssendProcEqns_
 
EqnBuffersendEqns_
 
std::vector< double > sendEqnSoln_
 
EqnBufferessBCEqns_
 
MPI_Comm comm_
 

Detailed Description

The EqnCommMgr (Equation communication manager) class is responsible for keeping track of equations that require communications. There are two types of equations in this class:

  1. Local equations which remote processors contribute to (e.g., because they share some local active nodes).
  2. Remote equations that the local processor contributes to (the mirror of case 1.).

Usage Notes:

  1. You can't call exchangeEqns until after exchangeIndices has been called.
  2. You can't call addSolnValues until after exchangeIndices has been called.

In general, usage will proceed like this:

in snl_fei::Structure::initComplete: addLocalEqn (a bunch of times, probably) addRemoteIndices (also a bunch of times)

exchangeIndices

getNumLocalEqns localEqnNumbers localIndicesPtr

in Filter::sumInElem/sumInElemMatrix/sumInElemRHS addRemoteEqn and/or addRemoteRHS

in Filter::exchangeRemoteEquations exchangeEqns getNumLocalEqns recvEqnNumbersPtr localIndicesPtr localCoefsPtr localRHSsPtr

in Filter::unpackSolution getNumLocalEqns localEqnNumbers addSolnValues exchangeSoln

in Filter various getSoln functions getNumSendEqns sendEqnNumbersPtr sendEqnSolnPtr

This class also provides general support functions for exchanging equations among processors. At the start of an all-to-all exchange, all processors usually know which equations they need to receive, and which processors they'll be receiving from, OR they know which equations they need to send, and which processors they'll be sending to. Usually they don't know both the sending and receiving information though. This eqn-processor pairing information (for either the "recv-equations" or the "send-equations") is held in a ProcEqns object. Thus, an all-to-all exchange of equation data requires two ProcEqns objects – one holding the recv info, the other holding the send info.

An EqnCommMgr function (mirrorProcEqns) is provided which does the communications necessary to populate the 'send' ProcEqns object given a populated 'recv' ProcEqns object, or vice-versa. Note that at this point, the equation-lengths need not be known by either the sending or the receiving processors.

Once the ProcEqns objects have been correctly mirrored, another EqnCommMgr function (mirrorProcEqnLengths) is available for mirroring the eqn-length data from one ProcEqns object to the other.

The next step is: given equations (with all associated data) that we need to send, along with previously known equation-numbers that we need to receive, an EqnCommMgr function (exchangeEqnBuffers) is provided to first send the equation-lengths to the receiving processors, followed by the actual equation data. At this point the exchange is complete. The equation data is supplied/returned in EqnBuffer objects.

Definition at line 104 of file fei_EqnCommMgr.hpp.

Constructor & Destructor Documentation

EqnCommMgr::EqnCommMgr ( MPI_Comm  comm,
bool  accumulate = true 
)

Constructor.

Parameters
localProcThe MPI rank of 'this' processor.

Definition at line 28 of file fei_EqnCommMgr.cpp.

References comm_, essBCEqns_, fei::localProc(), localProc_, recvEqns_, recvProcEqns_, sendEqns_, and sendProcEqns_.

Referenced by deepCopy().

EqnCommMgr::EqnCommMgr ( const EqnCommMgr src)

copy constructor

Definition at line 51 of file fei_EqnCommMgr.cpp.

EqnCommMgr::~EqnCommMgr ( )
virtual

Destructor.

Definition at line 68 of file fei_EqnCommMgr.cpp.

References essBCEqns_, recvEqns_, recvProcEqns_, sendEqns_, and sendProcEqns_.

Member Function Documentation

EqnCommMgr & EqnCommMgr::operator= ( const EqnCommMgr src)
EqnCommMgr * EqnCommMgr::deepCopy ( )
size_t EqnCommMgr::getNumSharingProcs ( )
inline

return the number of processors that share (contribute to) equations that are locally owned.

Definition at line 126 of file fei_EqnCommMgr.hpp.

References ProcEqns::getNumProcs(), and recvProcEqns_.

std::vector<int>& EqnCommMgr::sharingProcsPtr ( )
inline

Definition at line 127 of file fei_EqnCommMgr.hpp.

References ProcEqns::procsPtr(), and recvProcEqns_.

size_t EqnCommMgr::getNumOwnerProcs ( )
inline

return the number of processors that own equations that we share (equations that we contribute to).

Definition at line 132 of file fei_EqnCommMgr.hpp.

References ProcEqns::getNumProcs(), and sendProcEqns_.

std::vector<int>& EqnCommMgr::ownerProcsPtr ( )
inline

Definition at line 133 of file fei_EqnCommMgr.hpp.

References ProcEqns::procsPtr(), and sendProcEqns_.

void EqnCommMgr::addLocalEqn ( int  eqnNumber,
int  srcProc 
)

add a local equation to which a remote processor will be contributing.

Definition at line 131 of file fei_EqnCommMgr.cpp.

References ProcEqns::addEqn(), fei::console_out(), FEI_ENDL, localProc_, and recvProcEqns_.

Referenced by SNL_FEI_Structure::calculateSlaveEqns(), SNL_FEI_Structure::initializeEqnCommMgr(), and test_EqnCommMgr::test1().

void EqnCommMgr::addSolnValues ( int *  eqnNumbers,
double *  values,
int  num 
)
int EqnCommMgr::exchangeIndices ( std::ostream *  dbgOut = NULL)
int EqnCommMgr::exchangeEqns ( std::ostream *  dbgOut = NULL)
void EqnCommMgr::exchangeSoln ( )
int EqnCommMgr::mirrorProcEqns ( ProcEqns inProcEqns,
ProcEqns outProcEqns 
)

Support function to set up information needed for exchanging equation info among processors. Beginning assumption: we (the local processor) have a populated ProcEqns object (inProcEqns) which contains information pairing certain equations with certain remote processors. These can be equations that we will be receiving in an all-to-all exchange, or equations that we will be sending in an all-to-all exchange. In either case, the "mirror" of that information is needed before the exchange can be performed. i.e., if we know which eqns we'll be sending, then the mirror info concerns the eqns we'll be recv'ing, and vice-versa if we already know which eqns we'll be recv'ing.

This function is to obtain that mirror info, and return it in outProcEqns. Given a populated ProcEqns object, we want to populate the mirror ProcEqns object. Note that this function IGNORES any equation- lengths, if they are present. i.e., the eqn-length info in 'inProcEqns' is not referenced, and on completion the eqn-length info in 'outProcEqns' is not populated.

Parameters
inProcEqnsContains the input equation-numbers and associated processors.
outProcEqnsOutput.

Definition at line 744 of file fei_EqnCommMgr.cpp.

References ProcEqns::addEqn(), CHK_ERR, comm_, ProcEqns::eqnsPerProcPtr(), ERReturn, fei::exchangeData(), fei::mirrorProcs(), MPI_Request, MPI_SUCCESS, fei::numProcs(), ProcEqns::procEqnNumbersPtr(), and ProcEqns::procsPtr().

Referenced by exchangeIndices(), exchangeRemEssBCs(), gatherSharedBCs(), LinSysCoreFilter::getFromMatrix(), and LinSysCoreFilter::getFromRHS().

int EqnCommMgr::mirrorProcEqnLengths ( ProcEqns inProcEqns,
ProcEqns outProcEqns 
)

Support function to set up information needed for exchanging equation info among processors. This function plays a similar role to that of the above 'mirrorProcEqns' function, but exchanges the length information rather than the equation-number information. THIS FUNCTION ASSUMES that both inProcEqns and outProcEqns are already populated with eqn-number information.

Parameters
inProcEqnsContains equation-numbers, and their associated lengths.
outProcEqnsInput/Output. On entry, contains the equation-numbers (but not their lengths.
Returns
error non-zero if error occurs.

Definition at line 854 of file fei_EqnCommMgr.cpp.

References CHK_ERR, comm_, fei::exchangeData(), fei::numProcs(), ProcEqns::procEqnLengthsPtr(), and ProcEqns::procsPtr().

Referenced by exchangeIndices(), exchangeRemEssBCs(), gatherSharedBCs(), LinSysCoreFilter::getFromMatrix(), and LinSysCoreFilter::getFromRHS().

int EqnCommMgr::exchangeEqnBuffers ( MPI_Comm  comm,
ProcEqns sendProcEqns,
EqnBuffer sendEqns,
ProcEqns recvProcEqns,
EqnBuffer recvEqns,
bool  accumulate 
)
static
int EqnCommMgr::getNumLocalEqns ( )
inline
std::vector<int>& EqnCommMgr::localEqnNumbers ( )
inline
std::vector<fei::CSVec*>& EqnCommMgr::localEqns ( )
inline
std::vector<std::vector<double>*>* EqnCommMgr::localRHSsPtr ( )
inline

Definition at line 199 of file fei_EqnCommMgr.hpp.

References recvEqns_, and EqnBuffer::rhsCoefsPtr().

Referenced by LinSysCoreFilter::unpackRemoteContributions().

int EqnCommMgr::addRemoteEqn ( int  eqnNumber,
int  destProc,
const double *  coefs,
const int *  indices,
int  num 
)
int EqnCommMgr::addRemoteEqn ( int  eqnNumber,
const double *  coefs,
const int *  indices,
int  num 
)
int EqnCommMgr::addRemoteEqns ( fei::CSRMat mat,
bool  onlyIndices 
)
int EqnCommMgr::addRemoteRHS ( fei::CSVec vec,
int  rhsIndex 
)
void EqnCommMgr::setNumRHSs ( int  numRHSs)
int EqnCommMgr::addRemoteRHS ( int  eqnNumber,
int  destProc,
int  rhsIndex,
double  value 
)

Definition at line 907 of file fei_EqnCommMgr.cpp.

References EqnBuffer::addRHS(), EqnBuffer::newRHSData_, and sendEqns_.

int EqnCommMgr::addRemoteRHS ( int  eqnNumber,
int  rhsIndex,
double  value 
)

Definition at line 916 of file fei_EqnCommMgr.cpp.

References EqnBuffer::addRHS(), EqnBuffer::newRHSData_, and sendEqns_.

void EqnCommMgr::addRemoteIndices ( int  eqnNumber,
int  destProc,
int *  indices,
int  num 
)
int EqnCommMgr::getNumRemoteEqns ( )
inline

Definition at line 219 of file fei_EqnCommMgr.hpp.

References EqnBuffer::getNumEqns(), and sendEqns_.

std::vector<int>& EqnCommMgr::sendEqnNumbersPtr ( )
inline
double* EqnCommMgr::sendEqnSolnPtr ( )
inline
void EqnCommMgr::resetCoefs ( )
int EqnCommMgr::gatherSharedBCs ( EqnBuffer bcEqns)
int EqnCommMgr::exchangeRemEssBCs ( int *  essEqns,
int  numEssEqns,
double *  essAlpha,
double *  essGamma,
MPI_Comm  comm,
std::ostream *  dbgOut = NULL 
)
int EqnCommMgr::getNumRemEssBCEqns ( )
inline

Definition at line 233 of file fei_EqnCommMgr.hpp.

References essBCEqns_, and EqnBuffer::getNumEqns().

Referenced by LinSysCoreFilter::exchangeRemoteBCs().

std::vector<int>& EqnCommMgr::remEssBCEqnNumbersPtr ( )
inline

Definition at line 234 of file fei_EqnCommMgr.hpp.

References EqnBuffer::eqnNumbers(), and essBCEqns_.

Referenced by LinSysCoreFilter::exchangeRemoteBCs().

std::vector<fei::CSVec*>& EqnCommMgr::remEssBCEqns ( )
inline

Definition at line 235 of file fei_EqnCommMgr.hpp.

References EqnBuffer::eqns(), and essBCEqns_.

Referenced by LinSysCoreFilter::exchangeRemoteBCs().

int EqnCommMgr::exchangePtToBlkInfo ( snl_fei::PointBlockMap blkEqnMapper)
bool EqnCommMgr::newCoefData ( )
inline

Definition at line 239 of file fei_EqnCommMgr.hpp.

References EqnBuffer::newCoefData_, and recvEqns_.

Referenced by LinSysCoreFilter::unpackRemoteContributions().

bool EqnCommMgr::newRHSData ( )
inline

Definition at line 241 of file fei_EqnCommMgr.hpp.

References EqnBuffer::newRHSData_, and recvEqns_.

Referenced by LinSysCoreFilter::unpackRemoteContributions().

EqnBuffer* EqnCommMgr::getRecvEqns ( )
inline

Definition at line 246 of file fei_EqnCommMgr.hpp.

References recvEqns_.

Referenced by SNL_FEI_Structure::translateToReducedEqns().

EqnBuffer* EqnCommMgr::getSendEqns ( )
inline

Definition at line 247 of file fei_EqnCommMgr.hpp.

References sendEqns_.

Referenced by SNL_FEI_Structure::translateToReducedEqns().

ProcEqns* EqnCommMgr::getRecvProcEqns ( )
inline

Definition at line 248 of file fei_EqnCommMgr.hpp.

References recvProcEqns_.

Referenced by SNL_FEI_Structure::translateToReducedEqns().

ProcEqns* EqnCommMgr::getSendProcEqns ( )
inline

Definition at line 249 of file fei_EqnCommMgr.hpp.

References sendProcEqns_.

Referenced by SNL_FEI_Structure::translateToReducedEqns().

void EqnCommMgr::deleteEssBCs ( )
private
int EqnCommMgr::getSendProcNumber ( int  eqn)
private
int EqnCommMgr::consistencyCheck ( const char *  caller,
std::vector< int > &  recvProcs,
std::vector< int > &  recvProcTotalLengths,
std::vector< int > &  sendProcs,
std::vector< int > &  sendProcTotalLengths 
)
private

Member Data Documentation

bool EqnCommMgr::accumulate_
int EqnCommMgr::localProc_
private
ProcEqns* EqnCommMgr::recvProcEqns_
private
bool EqnCommMgr::exchangeIndicesCalled_
private

Definition at line 265 of file fei_EqnCommMgr.hpp.

Referenced by addSolnValues(), exchangeIndices(), and operator=().

EqnBuffer* EqnCommMgr::recvEqns_
private
std::vector<double> EqnCommMgr::solnValues_
private
ProcEqns* EqnCommMgr::sendProcEqns_
private
EqnBuffer* EqnCommMgr::sendEqns_
private
std::vector<double> EqnCommMgr::sendEqnSoln_
private

Definition at line 277 of file fei_EqnCommMgr.hpp.

Referenced by exchangeIndices(), exchangeSoln(), operator=(), and sendEqnSolnPtr().

EqnBuffer* EqnCommMgr::essBCEqns_
private
MPI_Comm EqnCommMgr::comm_
private

The documentation for this class was generated from the following files: