FEI Package Browser (Single Doxygen Collection)
Version of the Day
|
#include <fei_EqnCommMgr.hpp>
Public Member Functions | |
EqnCommMgr (MPI_Comm comm, bool accumulate=true) | |
EqnCommMgr (const EqnCommMgr &src) | |
EqnCommMgr & | operator= (const EqnCommMgr &src) |
virtual | ~EqnCommMgr () |
EqnCommMgr * | deepCopy () |
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 () |
EqnBuffer * | getRecvEqns () |
EqnBuffer * | getSendEqns () |
ProcEqns * | getRecvProcEqns () |
ProcEqns * | getSendProcEqns () |
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_ |
ProcEqns * | recvProcEqns_ |
bool | exchangeIndicesCalled_ |
EqnBuffer * | recvEqns_ |
std::vector< double > | solnValues_ |
ProcEqns * | sendProcEqns_ |
EqnBuffer * | sendEqns_ |
std::vector< double > | sendEqnSoln_ |
EqnBuffer * | essBCEqns_ |
MPI_Comm | comm_ |
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:
Usage Notes:
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.
EqnCommMgr::EqnCommMgr | ( | MPI_Comm | comm, |
bool | accumulate = true |
||
) |
Constructor.
localProc | The 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.
|
virtual |
Destructor.
Definition at line 68 of file fei_EqnCommMgr.cpp.
References essBCEqns_, recvEqns_, recvProcEqns_, sendEqns_, and sendProcEqns_.
EqnCommMgr & EqnCommMgr::operator= | ( | const EqnCommMgr & | src | ) |
assignment operator
Definition at line 79 of file fei_EqnCommMgr.cpp.
References accumulate_, EqnBuffer::deepCopy(), ProcEqns::deepCopy(), essBCEqns_, exchangeIndicesCalled_, recvEqns_, recvProcEqns_, sendEqns_, sendEqnSoln_, sendProcEqns_, and solnValues_.
EqnCommMgr * EqnCommMgr::deepCopy | ( | ) |
Produce a clone of 'this' object, including all of its internal data.
Definition at line 118 of file fei_EqnCommMgr.cpp.
References comm_, and EqnCommMgr().
Referenced by FEDataFilter::createEqnCommMgr_put(), LinSysCoreFilter::createEqnCommMgr_put(), FEDataFilter::FEDataFilter(), FEDataFilter::initialize(), LinSysCoreFilter::initialize(), LinSysCoreFilter::LinSysCoreFilter(), and test_EqnCommMgr::test1().
|
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_.
|
inline |
Definition at line 127 of file fei_EqnCommMgr.hpp.
References ProcEqns::procsPtr(), and recvProcEqns_.
|
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_.
|
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 | ||
) |
Definition at line 147 of file fei_EqnCommMgr.cpp.
References fei::binarySearch(), fei::console_out(), EqnBuffer::eqnNumbers(), exchangeIndicesCalled_, FEI_ENDL, recvEqns_, and solnValues_.
Referenced by test_EqnCommMgr::test1(), FEDataFilter::unpackSolution(), and LinSysCoreFilter::unpackSolution().
int EqnCommMgr::exchangeIndices | ( | std::ostream * | dbgOut = NULL | ) |
Definition at line 168 of file fei_EqnCommMgr.cpp.
References EqnBuffer::addIndices(), CHK_ERR, comm_, consistencyCheck(), fei::console_out(), dbgOut, EqnBuffer::eqnNumbers(), EqnBuffer::eqns(), ProcEqns::eqnsPerProcPtr(), ERReturn, exchangeIndicesCalled_, FEI_ENDL, FEI_OSTREAM, EqnBuffer::getEqnIndex(), EqnBuffer::getNumEqns(), ProcEqns::getNumProcs(), EqnBuffer::getNumRHSs(), mirrorProcEqnLengths(), mirrorProcEqns(), MPI_Request, MPI_SUCCESS, ProcEqns::procEqnLengthsPtr(), ProcEqns::procEqnNumbersPtr(), ProcEqns::procsPtr(), recvEqns_, recvProcEqns_, sendEqns_, sendEqnSoln_, sendProcEqns_, EqnBuffer::setNumRHSs(), ProcEqns::setProcEqnLengths(), fei::CSVec::size(), and solnValues_.
Referenced by SNL_FEI_Structure::formMatrixStructure(), SNL_FEI_Structure::initComplete(), and test_EqnCommMgr::test1().
int EqnCommMgr::exchangeEqns | ( | std::ostream * | dbgOut = NULL | ) |
Definition at line 438 of file fei_EqnCommMgr.cpp.
References accumulate_, CHK_ERR, comm_, dbgOut, exchangeEqnBuffers(), FEI_ENDL, FEI_OSTREAM, recvEqns_, recvProcEqns_, EqnBuffer::resetCoefs(), sendEqns_, and sendProcEqns_.
Referenced by LinSysCoreFilter::exchangeRemoteEquations().
void EqnCommMgr::exchangeSoln | ( | ) |
Definition at line 663 of file fei_EqnCommMgr.cpp.
References comm_, ProcEqns::eqnsPerProcPtr(), EqnBuffer::getEqnIndex(), EqnBuffer::getNumEqns(), ProcEqns::getNumProcs(), MPI_Comm, MPI_Request, ProcEqns::procEqnNumbersPtr(), ProcEqns::procsPtr(), recvEqns_, recvProcEqns_, sendEqns_, sendEqnSoln_, sendProcEqns_, and solnValues_.
Referenced by test_EqnCommMgr::test1(), FEDataFilter::unpackSolution(), and LinSysCoreFilter::unpackSolution().
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.
inProcEqns | Contains the input equation-numbers and associated processors. |
outProcEqns | Output. |
Definition at line 749 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().
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.
inProcEqns | Contains equation-numbers, and their associated lengths. |
outProcEqns | Input/Output. On entry, contains the equation-numbers (but not their lengths. |
Definition at line 859 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().
|
static |
Definition at line 466 of file fei_EqnCommMgr.cpp.
References EqnBuffer::addEqn(), EqnBuffer::addRHS(), CHK_ERR, EqnBuffer::eqns(), ProcEqns::eqnsPerProcPtr(), EqnBuffer::getEqnIndex(), ProcEqns::getNumProcs(), EqnBuffer::getNumRHSs(), MPI_Request, EqnBuffer::newCoefData_, EqnBuffer::newRHSData_, ProcEqns::procEqnLengthsPtr(), ProcEqns::procEqnNumbersPtr(), ProcEqns::procsPtr(), EqnBuffer::rhsCoefsPtr(), and EqnBuffer::setNumRHSs().
Referenced by exchangeEqns(), exchangeRemEssBCs(), gatherSharedBCs(), LinSysCoreFilter::getFromMatrix(), and LinSysCoreFilter::getFromRHS().
|
inline |
Definition at line 195 of file fei_EqnCommMgr.hpp.
References EqnBuffer::getNumEqns(), and recvEqns_.
Referenced by SNL_FEI_Structure::formMatrixStructure(), LinSysCoreFilter::unpackRemoteContributions(), FEDataFilter::unpackSolution(), and LinSysCoreFilter::unpackSolution().
|
inline |
Definition at line 197 of file fei_EqnCommMgr.hpp.
References EqnBuffer::eqnNumbers(), and recvEqns_.
Referenced by SNL_FEI_Structure::formMatrixStructure(), test_EqnCommMgr::test1(), LinSysCoreFilter::unpackRemoteContributions(), FEDataFilter::unpackSolution(), and LinSysCoreFilter::unpackSolution().
|
inline |
Definition at line 198 of file fei_EqnCommMgr.hpp.
References EqnBuffer::eqns(), and recvEqns_.
Referenced by SNL_FEI_Structure::formMatrixStructure(), and LinSysCoreFilter::unpackRemoteContributions().
|
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 | ||
) |
Definition at line 888 of file fei_EqnCommMgr.cpp.
References accumulate_, EqnBuffer::addEqn(), EqnBuffer::newCoefData_, and sendEqns_.
Referenced by addRemoteEqns(), LinSysCoreFilter::giveToBlkMatrix_symm_noSlaves(), LinSysCoreFilter::giveToMatrix(), LinSysCoreFilter::giveToMatrix_symm_noSlaves(), LinSysCoreFilter::storeNodalSendEqn(), and LinSysCoreFilter::storePenNodeSendData().
int EqnCommMgr::addRemoteEqn | ( | int | eqnNumber, |
const double * | coefs, | ||
const int * | indices, | ||
int | num | ||
) |
Definition at line 897 of file fei_EqnCommMgr.cpp.
References accumulate_, EqnBuffer::addEqn(), EqnBuffer::newCoefData_, and sendEqns_.
int EqnCommMgr::addRemoteEqns | ( | fei::CSRMat & | mat, |
bool | onlyIndices | ||
) |
Definition at line 1239 of file fei_EqnCommMgr.cpp.
References addRemoteEqn(), addRemoteIndices(), CHK_ERR, fei::CSRMat::getGraph(), fei::CSRMat::getPackedCoefs(), getSendProcNumber(), localProc_, fei::SparseRowGraph::packedColumnIndices, fei::SparseRowGraph::rowNumbers, and fei::SparseRowGraph::rowOffsets.
Referenced by SNL_FEI_Structure::assembleReducedStructure().
int EqnCommMgr::addRemoteRHS | ( | fei::CSVec & | vec, |
int | rhsIndex | ||
) |
Definition at line 1268 of file fei_EqnCommMgr.cpp.
References CHK_ERR, fei::CSVec::coefs(), getSendProcNumber(), fei::CSVec::indices(), and localProc_.
Referenced by LinSysCoreFilter::giveToRHS(), and LinSysCoreFilter::storePenNodeSendData().
void EqnCommMgr::setNumRHSs | ( | int | numRHSs | ) |
Definition at line 906 of file fei_EqnCommMgr.cpp.
References recvEqns_, sendEqns_, and EqnBuffer::setNumRHSs().
Referenced by SNL_FEI_Structure::calculateSlaveEqns(), SNL_FEI_Structure::initComplete(), FEDataFilter::initialize(), LinSysCoreFilter::initialize(), FEDataFilter::setNumRHSVectors(), LinSysCoreFilter::setNumRHSVectors(), SNL_FEI_Structure::SNL_FEI_Structure(), and test_EqnCommMgr::test1().
int EqnCommMgr::addRemoteRHS | ( | int | eqnNumber, |
int | destProc, | ||
int | rhsIndex, | ||
double | value | ||
) |
Definition at line 912 of file fei_EqnCommMgr.cpp.
References EqnBuffer::addRHS(), EqnBuffer::newRHSData_, and sendEqns_.
int EqnCommMgr::addRemoteRHS | ( | int | eqnNumber, |
int | rhsIndex, | ||
double | value | ||
) |
Definition at line 921 of file fei_EqnCommMgr.cpp.
References EqnBuffer::addRHS(), EqnBuffer::newRHSData_, and sendEqns_.
void EqnCommMgr::addRemoteIndices | ( | int | eqnNumber, |
int | destProc, | ||
int * | indices, | ||
int | num | ||
) |
Definition at line 928 of file fei_EqnCommMgr.cpp.
References ProcEqns::addEqn(), EqnBuffer::addIndices(), fei::console_out(), FEI_ENDL, sendEqns_, and sendProcEqns_.
Referenced by addRemoteEqns(), SNL_FEI_Structure::calculateSlaveEqns(), SNL_FEI_Structure::createBlkSymmEqnStructure(), SNL_FEI_Structure::createSymmEqnStructure(), SNL_FEI_Structure::initializeEqnCommMgr(), SNL_FEI_Structure::storeElementScatterBlkIndices_noSlaves(), SNL_FEI_Structure::storeElementScatterIndices(), SNL_FEI_Structure::storeElementScatterIndices_noSlaves(), SNL_FEI_Structure::storeNodalSendIndex(), SNL_FEI_Structure::storeNodalSendIndices(), and test_EqnCommMgr::test1().
|
inline |
Definition at line 219 of file fei_EqnCommMgr.hpp.
References EqnBuffer::getNumEqns(), and sendEqns_.
|
inline |
Definition at line 221 of file fei_EqnCommMgr.hpp.
References EqnBuffer::eqnNumbers(), and sendEqns_.
Referenced by FEDataFilter::getSharedRemoteSolnEntry(), and LinSysCoreFilter::getSharedRemoteSolnEntry().
|
inline |
Definition at line 223 of file fei_EqnCommMgr.hpp.
References sendEqnSoln_.
Referenced by FEDataFilter::getSharedRemoteSolnEntry(), and LinSysCoreFilter::getSharedRemoteSolnEntry().
void EqnCommMgr::resetCoefs | ( | ) |
Definition at line 942 of file fei_EqnCommMgr.cpp.
References essBCEqns_, EqnBuffer::getNumEqns(), EqnBuffer::newCoefData_, EqnBuffer::newRHSData_, recvEqns_, EqnBuffer::resetCoefs(), sendEqns_, and solnValues_.
Referenced by FEDataFilter::createEqnCommMgr_put(), LinSysCoreFilter::createEqnCommMgr_put(), LinSysCoreFilter::exchangeRemoteEquations(), FEDataFilter::resetMatrix(), LinSysCoreFilter::resetMatrix(), FEDataFilter::resetRHSVector(), LinSysCoreFilter::resetRHSVector(), and test_EqnCommMgr::test1().
int EqnCommMgr::gatherSharedBCs | ( | EqnBuffer & | bcEqns | ) |
Definition at line 959 of file fei_EqnCommMgr.cpp.
References ProcEqns::addEqn(), EqnBuffer::addEqn(), EqnBuffer::addEqns(), fei::binarySearch(), CHK_ERR, comm_, EqnBuffer::eqnNumbers(), EqnBuffer::eqns(), exchangeEqnBuffers(), mirrorProcEqnLengths(), mirrorProcEqns(), fei::numProcs(), ProcEqns::procEqnNumbersPtr(), ProcEqns::procsPtr(), sendEqns_, and sendProcEqns_.
Referenced by LinSysCoreFilter::implementAllBCs().
int EqnCommMgr::exchangeRemEssBCs | ( | int * | essEqns, |
int | numEssEqns, | ||
double * | essAlpha, | ||
double * | essGamma, | ||
MPI_Comm | comm, | ||
std::ostream * | dbgOut = NULL |
||
) |
Definition at line 1027 of file fei_EqnCommMgr.cpp.
References ProcEqns::addEqn(), EqnBuffer::addEqn(), fei::binarySearch(), CHK_ERR, dbgOut, EqnBuffer::eqnNumbers(), EqnBuffer::eqns(), essBCEqns_, exchangeEqnBuffers(), FEI_ENDL, FEI_OSTREAM, EqnBuffer::getNumEqns(), getSendProcNumber(), mirrorProcEqnLengths(), mirrorProcEqns(), sendEqns_, and ProcEqns::setProcEqnLengths().
Referenced by LinSysCoreFilter::exchangeRemoteBCs().
|
inline |
Definition at line 233 of file fei_EqnCommMgr.hpp.
References essBCEqns_, and EqnBuffer::getNumEqns().
Referenced by LinSysCoreFilter::exchangeRemoteBCs().
|
inline |
Definition at line 234 of file fei_EqnCommMgr.hpp.
References EqnBuffer::eqnNumbers(), and essBCEqns_.
Referenced by LinSysCoreFilter::exchangeRemoteBCs().
|
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 | ) |
Definition at line 1146 of file fei_EqnCommMgr.cpp.
References CHK_ERR, comm_, EqnBuffer::eqns(), snl_fei::PointBlockMap::eqnToBlkEqn(), fei::exchangeData(), snl_fei::PointBlockMap::getBlkEqnSize(), ProcEqns::getNumProcs(), snl_fei::PointBlockMap::getPtEqns(), ProcEqns::procsPtr(), recvEqns_, recvProcEqns_, sendEqns_, sendProcEqns_, and snl_fei::PointBlockMap::setEqn().
Referenced by SNL_FEI_Structure::initComplete().
|
inline |
Definition at line 239 of file fei_EqnCommMgr.hpp.
References EqnBuffer::newCoefData_, and recvEqns_.
Referenced by LinSysCoreFilter::unpackRemoteContributions().
|
inline |
Definition at line 241 of file fei_EqnCommMgr.hpp.
References EqnBuffer::newRHSData_, and recvEqns_.
Referenced by LinSysCoreFilter::unpackRemoteContributions().
|
inline |
Definition at line 246 of file fei_EqnCommMgr.hpp.
References recvEqns_.
Referenced by SNL_FEI_Structure::translateToReducedEqns().
|
inline |
Definition at line 247 of file fei_EqnCommMgr.hpp.
References sendEqns_.
Referenced by SNL_FEI_Structure::translateToReducedEqns().
|
inline |
Definition at line 248 of file fei_EqnCommMgr.hpp.
References recvProcEqns_.
Referenced by SNL_FEI_Structure::translateToReducedEqns().
|
inline |
Definition at line 249 of file fei_EqnCommMgr.hpp.
References sendProcEqns_.
Referenced by SNL_FEI_Structure::translateToReducedEqns().
|
private |
|
private |
Definition at line 1285 of file fei_EqnCommMgr.cpp.
References ProcEqns::getNumProcs(), ProcEqns::procEqnNumbersPtr(), ProcEqns::procsPtr(), and sendProcEqns_.
Referenced by addRemoteEqns(), addRemoteRHS(), and exchangeRemEssBCs().
|
private |
Definition at line 363 of file fei_EqnCommMgr.cpp.
References fei::Allgatherv(), CHK_ERR, comm_, fei::console_out(), FEI_ENDL, fei::GlobalSum(), and localProc_.
Referenced by exchangeIndices().
bool EqnCommMgr::accumulate_ |
Definition at line 244 of file fei_EqnCommMgr.hpp.
Referenced by addRemoteEqn(), FEDataFilter::createEqnCommMgr_put(), LinSysCoreFilter::createEqnCommMgr_put(), exchangeEqns(), and operator=().
|
private |
Definition at line 261 of file fei_EqnCommMgr.hpp.
Referenced by addLocalEqn(), addRemoteEqns(), addRemoteRHS(), consistencyCheck(), and EqnCommMgr().
|
private |
Definition at line 263 of file fei_EqnCommMgr.hpp.
Referenced by addLocalEqn(), EqnCommMgr(), exchangeEqns(), exchangeIndices(), exchangePtToBlkInfo(), exchangeSoln(), getNumSharingProcs(), getRecvProcEqns(), operator=(), sharingProcsPtr(), and ~EqnCommMgr().
|
private |
Definition at line 265 of file fei_EqnCommMgr.hpp.
Referenced by addSolnValues(), exchangeIndices(), and operator=().
|
private |
Definition at line 268 of file fei_EqnCommMgr.hpp.
Referenced by addSolnValues(), EqnCommMgr(), exchangeEqns(), exchangeIndices(), exchangePtToBlkInfo(), exchangeSoln(), getNumLocalEqns(), getRecvEqns(), localEqnNumbers(), localEqns(), localRHSsPtr(), newCoefData(), newRHSData(), operator=(), resetCoefs(), setNumRHSs(), and ~EqnCommMgr().
|
private |
Definition at line 270 of file fei_EqnCommMgr.hpp.
Referenced by addSolnValues(), exchangeIndices(), exchangeSoln(), operator=(), and resetCoefs().
|
private |
Definition at line 273 of file fei_EqnCommMgr.hpp.
Referenced by addRemoteIndices(), EqnCommMgr(), exchangeEqns(), exchangeIndices(), exchangePtToBlkInfo(), exchangeSoln(), gatherSharedBCs(), getNumOwnerProcs(), getSendProcEqns(), getSendProcNumber(), operator=(), ownerProcsPtr(), and ~EqnCommMgr().
|
private |
Definition at line 275 of file fei_EqnCommMgr.hpp.
Referenced by addRemoteEqn(), addRemoteIndices(), addRemoteRHS(), EqnCommMgr(), exchangeEqns(), exchangeIndices(), exchangePtToBlkInfo(), exchangeRemEssBCs(), exchangeSoln(), gatherSharedBCs(), getNumRemoteEqns(), getSendEqns(), operator=(), resetCoefs(), sendEqnNumbersPtr(), setNumRHSs(), and ~EqnCommMgr().
|
private |
Definition at line 277 of file fei_EqnCommMgr.hpp.
Referenced by exchangeIndices(), exchangeSoln(), operator=(), and sendEqnSolnPtr().
|
private |
Definition at line 282 of file fei_EqnCommMgr.hpp.
Referenced by EqnCommMgr(), exchangeRemEssBCs(), getNumRemEssBCEqns(), operator=(), remEssBCEqnNumbersPtr(), remEssBCEqns(), resetCoefs(), and ~EqnCommMgr().
|
private |
Definition at line 284 of file fei_EqnCommMgr.hpp.
Referenced by consistencyCheck(), deepCopy(), EqnCommMgr(), exchangeEqns(), exchangeIndices(), exchangePtToBlkInfo(), exchangeSoln(), gatherSharedBCs(), mirrorProcEqnLengths(), and mirrorProcEqns().