Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT > Class Template Reference

#include <Panzer_BlockedEpetraLinearObjFactory.hpp>

Inheritance diagram for panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >:
Inheritance graph
[legend]

Classes

class  DOFManagerContainer
 

Public Member Functions

 BlockedEpetraLinearObjFactory (const Teuchos::RCP< const Teuchos::MpiComm< int > > &comm, const Teuchos::RCP< const GlobalIndexer > &gidProvider, bool useDiscreteAdjoint=false)
 
 BlockedEpetraLinearObjFactory (const Teuchos::RCP< const Teuchos::MpiComm< int > > &comm, const Teuchos::RCP< const GlobalIndexer > &gidProvider, const Teuchos::RCP< const GlobalIndexer > &colGidProvider, bool useDiscreteAdjoint=false)
 
virtual ~BlockedEpetraLinearObjFactory ()
 
virtual void readVector (const std::string &identifier, LinearObjContainer &loc, int id) const
 
virtual void writeVector (const std::string &identifier, const LinearObjContainer &loc, int id) const
 
virtual Teuchos::RCP
< LinearObjContainer
buildLinearObjContainer () const
 
virtual Teuchos::RCP
< LinearObjContainer
buildPrimitiveLinearObjContainer () const
 
virtual Teuchos::RCP
< LinearObjContainer
buildGhostedLinearObjContainer () const
 
virtual Teuchos::RCP
< LinearObjContainer
buildPrimitiveGhostedLinearObjContainer () const
 
virtual void globalToGhostContainer (const LinearObjContainer &container, LinearObjContainer &ghostContainer, int) const
 
virtual void ghostToGlobalContainer (const LinearObjContainer &ghostContainer, LinearObjContainer &container, int) const
 
virtual void adjustForDirichletConditions (const LinearObjContainer &localBCRows, const LinearObjContainer &globalBCRows, LinearObjContainer &ghostedObjs, bool zeroVectorRows=false, bool adjustX=false) const
 
virtual void applyDirichletBCs (const LinearObjContainer &counter, LinearObjContainer &result) const
 
virtual Teuchos::RCP
< ReadOnlyVector_GlobalEvaluationData
buildReadOnlyDomainContainer () const
 
virtual Teuchos::RCP
< WriteVector_GlobalEvaluationData
buildWriteDomainContainer () const
 
virtual Teuchos::MpiComm< int > getComm () const
 
template<typename EvalT >
Teuchos::RCP
< panzer::CloneableEvaluator
buildScatter () const
 Use preconstructed scatter evaluators. More...
 
template<typename EvalT >
Teuchos::RCP
< panzer::CloneableEvaluator
buildGather () const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP
< panzer::CloneableEvaluator
buildGatherTangent () const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP
< panzer::CloneableEvaluator
buildGatherDomain () const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP
< panzer::CloneableEvaluator
buildGatherOrientation () const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP
< panzer::CloneableEvaluator
buildScatterDirichlet () const
 Use preconstructed dirichlet scatter evaluators. More...
 
void initializeContainer (int, LinearObjContainer &loc) const
 
void initializeGhostedContainer (int, LinearObjContainer &loc) const
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
getThyraDomainSpace () const
 Get the domain vector space (x and dxdt) More...
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
getThyraRangeSpace () const
 Get the range vector space (f) More...
 
Teuchos::RCP
< Thyra::VectorBase< double > > 
getThyraDomainVector () const
 Get a domain vector. More...
 
Teuchos::RCP
< Thyra::VectorBase< double > > 
getThyraRangeVector () const
 Get a range vector. More...
 
Teuchos::RCP
< Thyra::LinearOpBase< double > > 
getThyraMatrix () const
 Get a Thyra operator. More...
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
getGhostedThyraDomainSpace () const
 Get the domain vector space (x and dxdt) More...
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
getGhostedThyraDomainSpace2 () const
 Get or create the ghosted Thyra domain space. More...
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
getGhostedThyraRangeSpace () const
 Get the range vector space (f) More...
 
Teuchos::RCP
< Thyra::VectorBase< double > > 
getGhostedThyraDomainVector () const
 Get a domain vector. More...
 
Teuchos::RCP
< Thyra::VectorBase< double > > 
getGhostedThyraRangeVector () const
 Get a range vector. More...
 
Teuchos::RCP
< Thyra::LinearOpBase< double > > 
getGhostedThyraMatrix () const
 Get a Thyra operator. More...
 
virtual const Teuchos::RCP
< Epetra_Map
getMap (int i) const
 get the map from the matrix More...
 
virtual const Teuchos::RCP
< Epetra_Map
getColMap (int i) const
 get the map from the matrix More...
 
virtual const Teuchos::RCP
< Epetra_Map
getGhostedMap (int i) const
 get the ghosted map from the matrix More...
 
virtual const Teuchos::RCP
< Epetra_Map
getGhostedMap2 (int i) const
 Get or create the i-th ghosted map. More...
 
virtual const Teuchos::RCP
< Epetra_Map
getGhostedColMap (int i) const
 get the ghosted map from the matrix More...
 
virtual const Teuchos::RCP
< Epetra_Map
getGhostedColMap2 (int i) const
 Get or create the i-th ghosted column map. More...
 
virtual const Teuchos::RCP
< Epetra_CrsGraph
getGraph (int i, int j) const
 get the graph of the crs matrix More...
 
virtual const Teuchos::RCP
< Epetra_CrsGraph
getGhostedGraph (int i, int j) const
 get the ghosted graph of the crs matrix More...
 
virtual const Teuchos::RCP
< Epetra_Import
getGhostedImport (int i) const
 get importer for converting an overalapped object to a "normal" object More...
 
virtual const Teuchos::RCP
< Epetra_Import
getGhostedImport2 (int i) const
 Get or create the i-th ghosted importer corresponding to the i-th ghosted map. More...
 
virtual const Teuchos::RCP
< Epetra_Import
getGhostedColImport (int i) const
 get importer for converting an overalapped object to a "normal" object More...
 
virtual const Teuchos::RCP
< Epetra_Import
getGhostedColImport2 (int i) const
 Get or create the i-th ghosted column importer corresponding to the i-th ghosted column map. More...
 
virtual const Teuchos::RCP
< Epetra_Export
getGhostedExport (int j) const
 get exporter for converting an overalapped object to a "normal" object More...
 
virtual const Teuchos::RCP
< Epetra_Export
getGhostedExport2 (int i) const
 Get or create the i-th ghosted exporter corresponding to the i-th ghosted map. More...
 
virtual const Teuchos::RCP
< Epetra_Export
getGhostedColExport (int j) const
 get exporter for converting an overalapped object to a "normal" object More...
 
virtual const Teuchos::RCP
< Epetra_Export
getGhostedColExport2 (int i) const
 Get or create the i-th ghosted column exporter corresponding to the i-th ghosted column map. More...
 
virtual const Teuchos::RCP
< const Epetra_Comm
getEpetraComm () const
 get exporter for converting an overalapped object to a "normal" object More...
 
Teuchos::RCP< Epetra_CrsMatrixgetEpetraMatrix (int i, int j) const
 
Teuchos::RCP< Epetra_CrsMatrixgetGhostedEpetraMatrix (int i, int j) const
 
int getBlockRowCount () const
 how many block rows More...
 
int getBlockColCount () const
 how many block columns More...
 
Teuchos::RCP< const
panzer::BlockedDOFManager
getGlobalIndexer () const
 
Teuchos::RCP< const
panzer::GlobalIndexer
getRangeGlobalIndexer () const
 Get the range global indexer object associated with this factory. More...
 
Teuchos::RCP< const
panzer::GlobalIndexer
getDomainGlobalIndexer () const
 Get the domain global indexer object associated with this factory. More...
 
const std::vector
< Teuchos::RCP< const
GlobalIndexer > > & 
getRangeGlobalIndexers () const
 Get global indexers associated with the blocks. More...
 
const std::vector
< Teuchos::RCP< const
GlobalIndexer > > & 
getDomainGlobalIndexers () const
 Get global indexers associated with the blocks. More...
 
void addExcludedPair (int rowBlock, int colBlock)
 exclude a block pair from the matrix More...
 
void addExcludedPairs (const std::vector< std::pair< int, int > > &exPairs)
 exclude a vector of pairs from the matrix More...
 
- Public Member Functions inherited from panzer::LinearObjFactory< Traits >
virtual ~LinearObjFactory ()
 
template<typename BuilderT >
void buildGatherScatterEvaluators (const BuilderT &builder)
 
template<typename EvalT >
Teuchos::RCP< PHX::Evaluator
< Traits > > 
buildScatter (const Teuchos::ParameterList &pl) const
 Use preconstructed scatter evaluators. More...
 
template<typename EvalT >
Teuchos::RCP< PHX::Evaluator
< Traits > > 
buildGather (const Teuchos::ParameterList &pl) const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP< PHX::Evaluator
< Traits > > 
buildGatherTangent (const Teuchos::ParameterList &pl) const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP< PHX::Evaluator
< Traits > > 
buildGatherDomain (const Teuchos::ParameterList &pl) const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP< PHX::Evaluator
< Traits > > 
buildGatherOrientation (const Teuchos::ParameterList &pl) const
 Use preconstructed gather evaluators. More...
 
template<typename EvalT >
Teuchos::RCP< PHX::Evaluator
< Traits > > 
buildScatterDirichlet (const Teuchos::ParameterList &pl) const
 Use preconstructed dirichlet scatter evaluators. More...
 
virtual void beginFill (LinearObjContainer &) const
 
virtual void endFill (LinearObjContainer &) const
 
- Public Member Functions inherited from panzer::ThyraObjFactory< double >
virtual ~ThyraObjFactory ()
 

Protected Member Functions

void initializeContainer_internal (int mem, ThyraObjContainer< double > &loc) const
 
void initializeGhostedContainer_internal (int mem, ThyraObjContainer< double > &loc) const
 
Teuchos::RCP< const GlobalIndexergetGlobalIndexer (int i) const
 
Teuchos::RCP< const GlobalIndexergetColGlobalIndexer (int i) const
 
void makeRoomForBlocks (std::size_t blockCnt, std::size_t colBlockCnt=0)
 Allocate the space in the std::vector objects so we can fill with appropriate Epetra data. More...
 
void ghostToGlobalThyraVector (const Teuchos::RCP< const Thyra::VectorBase< double > > &in, const Teuchos::RCP< Thyra::VectorBase< double > > &out, bool col) const
 
void ghostToGlobalThyraMatrix (const Thyra::LinearOpBase< double > &in, Thyra::LinearOpBase< double > &out) const
 
void globalToGhostThyraVector (const Teuchos::RCP< const Thyra::VectorBase< double > > &in, const Teuchos::RCP< Thyra::VectorBase< double > > &out, bool col) const
 
void adjustForDirichletConditions (const Epetra_Vector &local_bcs, const Epetra_Vector &global_bcs, const Teuchos::Ptr< Epetra_Vector > &f, const Teuchos::Ptr< Epetra_CrsMatrix > &A, bool zeroVectorRows) const
 
void ghostToGlobalEpetraVector (int i, const Epetra_Vector &in, Epetra_Vector &out, bool col) const
 
void globalToGhostEpetraVector (int i, const Epetra_Vector &in, Epetra_Vector &out, bool col) const
 
void ghostToGlobalEpetraMatrix (int blockRow, const Epetra_CrsMatrix &in, Epetra_CrsMatrix &out) const
 
virtual const Teuchos::RCP
< Epetra_Map
buildMap (int i) const
 Build the i-th owned map from the owned indices of the i-th global indexer. More...
 
virtual const Teuchos::RCP
< Epetra_Map
buildGhostedMap (int i) const
 Build the i-th ghosted map from the owned and ghosted indices of the i-th global indexer. More...
 
virtual const Teuchos::RCP
< Epetra_Map
buildGhostedMap2 (int i) const
 Build the i-th ghosted map from the ghosted indices of the i-th global indexer. More...
 
virtual const Teuchos::RCP
< Epetra_Map
buildColMap (int i) const
 Build the i-th owned column map from the owned indices of the i-th (column) global indexer. More...
 
virtual const Teuchos::RCP
< Epetra_Map
buildColGhostedMap (int i) const
 Build the i-th ghosted column map from the owned and ghosted indices of the i-th (column) global indexer. More...
 
virtual const Teuchos::RCP
< Epetra_Map
buildColGhostedMap2 (int i) const
 Build the i-th ghosted column map from the ghosted indices of the i-th (column) global indexer. More...
 
virtual const Teuchos::RCP
< Epetra_CrsGraph
buildGraph (int i, int j) const
 
virtual const Teuchos::RCP
< Epetra_CrsGraph
buildGhostedGraph (int i, int j, bool optimizeStorage) const
 
virtual const Teuchos::RCP
< Epetra_CrsGraph
buildFilteredGhostedGraph (int i, int j) const
 

Protected Attributes

Teuchos::RCP< const
DOFManagerContainer
rowDOFManagerContainer_
 
Teuchos::RCP< const
DOFManagerContainer
colDOFManagerContainer_
 
bool useColGidProviders_
 
std::unordered_set< std::pair
< int, int >
, panzer::pair_hash
excludedPairs_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
rangeSpace_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
domainSpace_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
ghostedRangeSpace_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< double > > 
ghostedDomainSpace_
 
Teuchos::RCP< const Epetra_CommeComm_
 
Teuchos::RCP< const
Teuchos::OpaqueWrapper
< MPI_Comm > > 
rawMpiComm_
 
Teuchos::RCP< Teuchos::MpiComm
< int > > 
tComm_
 
std::vector< Teuchos::RCP
< Epetra_Map > > 
maps_
 The list of owned maps corresponding to the owned indices of the global indexers. More...
 
std::vector< Teuchos::RCP
< Epetra_Map > > 
ghostedMaps_
 The list of ghosted maps corresponding to the owned and ghosted indices of the global indexers. More...
 
std::vector< Teuchos::RCP
< Epetra_Map > > 
ghostedMaps2_
 The list of ghosted maps corresponding to the ghosted indices of the global indexers. More...
 
std::vector< Teuchos::RCP
< Epetra_Import > > 
importers_
 The list of ghosted importers corresponding to ghostedMaps_. More...
 
std::vector< Teuchos::RCP
< Epetra_Import > > 
importers2_
 The list of ghosted importers corresponding to ghostedMaps2_. More...
 
std::vector< Teuchos::RCP
< Epetra_Export > > 
exporters_
 
std::vector< Teuchos::RCP
< Epetra_Map > > 
colMaps_
 The list of owned column maps corresponding to the owned indices of the (column) global indexers. More...
 
std::vector< Teuchos::RCP
< Epetra_Map > > 
colGhostedMaps_
 The list of ghosted column maps corresponding to the owned and ghosted indices of the (column) global indexers. More...
 
std::vector< Teuchos::RCP
< Epetra_Map > > 
colGhostedMaps2_
 The list of ghosted column maps corresponding to the ghosted indices of the (column) global indexers. More...
 
std::vector< Teuchos::RCP
< Epetra_Import > > 
colImporters_
 The list of ghosted importers corresponding to colGhostedMaps_. More...
 
std::vector< Teuchos::RCP
< Epetra_Import > > 
colImporters2_
 The list of ghosted importers corresponding to colGhostedMaps2_. More...
 
std::vector< Teuchos::RCP
< Epetra_Export > > 
colExporters_
 
std::unordered_map< std::pair
< int, int >, Teuchos::RCP
< Epetra_CrsGraph >
, panzer::pair_hash
graphs_
 
std::unordered_map< std::pair
< int, int >, Teuchos::RCP
< Epetra_CrsGraph >
, panzer::pair_hash
ghostedGraphs_
 
bool useDiscreteAdjoint_
 

Detailed Description

template<typename Traits, typename LocalOrdinalT>
class panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >

Definition at line 86 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

Constructor & Destructor Documentation

template<typename Traits , typename LocalOrdinalT >
panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::BlockedEpetraLinearObjFactory ( const Teuchos::RCP< const Teuchos::MpiComm< int > > &  comm,
const Teuchos::RCP< const GlobalIndexer > &  gidProvider,
bool  useDiscreteAdjoint = false 
)
template<typename Traits , typename LocalOrdinalT >
panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::BlockedEpetraLinearObjFactory ( const Teuchos::RCP< const Teuchos::MpiComm< int > > &  comm,
const Teuchos::RCP< const GlobalIndexer > &  gidProvider,
const Teuchos::RCP< const GlobalIndexer > &  colGidProvider,
bool  useDiscreteAdjoint = false 
)
template<typename Traits , typename LocalOrdinalT >
panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::~BlockedEpetraLinearObjFactory ( )
virtual

Member Function Documentation

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::readVector ( const std::string &  identifier,
LinearObjContainer loc,
int  id 
) const
virtual

Read in a vector from a file. Fill a particular vector in the linear object container.

Parameters
[in]identifierKey for specifying which file(s) to read
[in]locLinear object container to fill with the vector
[in]idId for the field to be filled

Implements panzer::LinearObjFactory< Traits >.

Definition at line 143 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::writeVector ( const std::string &  identifier,
const LinearObjContainer loc,
int  id 
) const
virtual

Write in a vector from a file. Fill a particular vector in the linear object container.

Parameters
[in]identifierKey for specifying which file(s) to read
[in]locLinear object container to fill with the vector
[in]idId for the field to be filled

Implements panzer::LinearObjFactory< Traits >.

Definition at line 193 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< LinearObjContainer > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildLinearObjContainer ( ) const
virtual

Build a container with all the neccessary linear algebra objects. This is the non-ghosted version.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 237 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits, typename LocalOrdinalT>
virtual Teuchos::RCP<LinearObjContainer> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildPrimitiveLinearObjContainer ( ) const
inlinevirtual

Build a container with all the neccessary linear algebra objects, purely on the single physics. This gives linear algebra objects that are relevant for a single physics solve. In many cases this is simply a call to buildLinearObjContainer however, in a few important cases (for instance in stochastic galerkin methods) this will return a container for a single instantiation of the physics. This is the non-ghosted version.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 109 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< LinearObjContainer > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGhostedLinearObjContainer ( ) const
virtual

Build a container with all the neccessary linear algebra objects. This is the ghosted version.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 258 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits, typename LocalOrdinalT>
virtual Teuchos::RCP<LinearObjContainer> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildPrimitiveGhostedLinearObjContainer ( ) const
inlinevirtual

Build a container with all the neccessary linear algebra objects, purely on the single physics. This gives linear algebra objects that are relevant for a single physics solve. In many cases this is simply a call to buildGhostedLinearObjContainer however, in a few important cases (for instance in stochastic galerkin methods) this will return a container for a single instantiation of the physics. This is the ghosted version.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 114 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::globalToGhostContainer ( const LinearObjContainer container,
LinearObjContainer ghostContainer,
int  mem 
) const
virtual
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostToGlobalContainer ( const LinearObjContainer ghostContainer,
LinearObjContainer container,
int  mem 
) const
virtual
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::adjustForDirichletConditions ( const LinearObjContainer localBCRows,
const LinearObjContainer globalBCRows,
LinearObjContainer ghostedObjs,
bool  zeroVectorRows = false,
bool  adjustX = false 
) const
virtual

Adjust the residual vector and Jacobian matrix (if they exist) for applied dirichlet conditions. The adjustment considers if a boundary condition was set globally and locally and based on that result adjust the ghosted matrix and residual vector so that when they are summed across processors they resulting Dirichlet condition is correct.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 364 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::applyDirichletBCs ( const LinearObjContainer counter,
LinearObjContainer result 
) const
virtual

Adjust a vector by replacing selected rows with the value of the evaluated dirichlet conditions. This is handled through the standard container mechanism.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 501 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< ReadOnlyVector_GlobalEvaluationData > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildReadOnlyDomainContainer ( ) const
virtual

Build a GlobalEvaluationDataContainer that handles all domain communication. This is used primarily for gather operations and hides the allocation and usage of the ghosted vector from the user.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 545 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< WriteVector_GlobalEvaluationData > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildWriteDomainContainer ( ) const
virtual

Build a GlobalEvaluationDataContainer that handles all domain communication. This is used primarily for gather operations and hides the allocation and usage of the ghosted vector from the user.

Definition at line 587 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::MpiComm< int > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getComm ( ) const
virtual

Acess to the MPI Comm used in constructing this LOF.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 623 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits, typename LocalOrdinalT>
template<typename EvalT >
Teuchos::RCP<panzer::CloneableEvaluator> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildScatter ( ) const
inline

Use preconstructed scatter evaluators.

Definition at line 155 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
template<typename EvalT >
Teuchos::RCP<panzer::CloneableEvaluator > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGather ( ) const
inline

Use preconstructed gather evaluators.

Definition at line 170 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
template<typename EvalT >
Teuchos::RCP<panzer::CloneableEvaluator > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGatherTangent ( ) const
inline

Use preconstructed gather evaluators.

Definition at line 180 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
template<typename EvalT >
Teuchos::RCP<panzer::CloneableEvaluator > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGatherDomain ( ) const
inline

Use preconstructed gather evaluators.

Definition at line 190 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
template<typename EvalT >
Teuchos::RCP<panzer::CloneableEvaluator > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGatherOrientation ( ) const
inline

Use preconstructed gather evaluators.

Definition at line 199 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
template<typename EvalT >
Teuchos::RCP<panzer::CloneableEvaluator> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildScatterDirichlet ( ) const
inline

Use preconstructed dirichlet scatter evaluators.

Definition at line 204 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::initializeContainer ( int  mem,
LinearObjContainer loc 
) const
virtual

Initialize container with a specific set of member values.

Note
This will overwrite everything in the container and zero out values not requested.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 630 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::initializeGhostedContainer ( int  mem,
LinearObjContainer loc 
) const
virtual

Initialize container with a specific set of member values.

Note
This will overwrite everything in the container and zero out values not requested.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 640 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getThyraDomainSpace ( ) const
virtual

Get the domain vector space (x and dxdt)

Implements panzer::ThyraObjFactory< double >.

Definition at line 778 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getThyraRangeSpace ( ) const
virtual

Get the range vector space (f)

Implements panzer::ThyraObjFactory< double >.

Definition at line 801 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Thyra::VectorBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getThyraDomainVector ( ) const

Get a domain vector.

Definition at line 824 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Thyra::VectorBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getThyraRangeVector ( ) const

Get a range vector.

Definition at line 835 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Thyra::LinearOpBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getThyraMatrix ( ) const
virtual

Get a Thyra operator.

Implements panzer::ThyraObjFactory< double >.

Definition at line 846 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedThyraDomainSpace ( ) const

Get the domain vector space (x and dxdt)

Definition at line 888 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedThyraDomainSpace2 ( ) const

Get or create the ghosted Thyra domain space.

Get the vector space corresponding to the ghosted domain. If it does not yet exist, create it from the ghosted column map(s).

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Returns
The vector space corresponding to the ghosted domain.

Definition at line 923 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Thyra::VectorSpaceBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedThyraRangeSpace ( ) const

Get the range vector space (f)

Definition at line 952 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Thyra::VectorBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedThyraDomainVector ( ) const

Get a domain vector.

Definition at line 975 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Thyra::VectorBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedThyraRangeVector ( ) const

Get a range vector.

Definition at line 986 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Thyra::LinearOpBase< double > > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedThyraMatrix ( ) const

Get a Thyra operator.

Definition at line 997 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getMap ( int  i) const
virtual

get the map from the matrix

Definition at line 1197 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getColMap ( int  i) const
virtual

get the map from the matrix

Definition at line 1208 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedMap ( int  i) const
virtual

get the ghosted map from the matrix

Definition at line 1227 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedMap2 ( int  i) const
virtual

Get or create the i-th ghosted map.

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Parameters
[in]iThe index into the list of ghosted maps.
Returns
The i-th ghosted map.

Definition at line 1243 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedColMap ( int  i) const
virtual

get the ghosted map from the matrix

Definition at line 1259 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedColMap2 ( int  i) const
virtual

Get or create the i-th ghosted column map.

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Parameters
[in]iThe index into the list of ghosted column maps.
Returns
The i-th ghosted column map.

Definition at line 1277 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_CrsGraph > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGraph ( int  i,
int  j 
) const
virtual

get the graph of the crs matrix

Definition at line 1290 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_CrsGraph > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedGraph ( int  i,
int  j 
) const
virtual

get the ghosted graph of the crs matrix

Definition at line 1309 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Import > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedImport ( int  i) const
virtual

get importer for converting an overalapped object to a "normal" object

Definition at line 1334 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Import > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedImport2 ( int  i) const
virtual

Get or create the i-th ghosted importer corresponding to the i-th ghosted map.

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Parameters
[in]iThe index into the list of ghosted importers.
Returns
The i-th ghosted importer.

Definition at line 1351 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Import > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedColImport ( int  i) const
virtual

get importer for converting an overalapped object to a "normal" object

Definition at line 1368 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Import > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedColImport2 ( int  i) const
virtual

Get or create the i-th ghosted column importer corresponding to the i-th ghosted column map.

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Parameters
[in]iThe index into the list of ghosted column importers.
Returns
The i-th ghosted column importer.

Definition at line 1388 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Export > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedExport ( int  j) const
virtual

get exporter for converting an overalapped object to a "normal" object

Definition at line 1408 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Export > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedExport2 ( int  i) const
virtual

Get or create the i-th ghosted exporter corresponding to the i-th ghosted map.

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Parameters
[in]iThe index into the list of ghosted exporters.
Returns
The i-th ghosted exporter.

Definition at line 1425 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Export > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedColExport ( int  j) const
virtual

get exporter for converting an overalapped object to a "normal" object

Definition at line 1442 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Export > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedColExport2 ( int  i) const
virtual

Get or create the i-th ghosted column exporter corresponding to the i-th ghosted column map.

Note
This "version 2" routine works with non-overlapping owned and ghosted maps.
Parameters
[in]iThe index into the list of ghosted column exporters.
Returns
The i-th ghosted column exporter.

Definition at line 1462 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< const Epetra_Comm > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getEpetraComm ( ) const
virtual

get exporter for converting an overalapped object to a "normal" object

Definition at line 1743 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Epetra_CrsMatrix > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getEpetraMatrix ( int  i,
int  j 
) const
template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Epetra_CrsMatrix > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGhostedEpetraMatrix ( int  i,
int  j 
) const
template<typename Traits , typename LocalOrdinalT >
int panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getBlockRowCount ( ) const

how many block rows

Definition at line 1750 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
int panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getBlockColCount ( ) const

how many block columns

Definition at line 1757 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const panzer::BlockedDOFManager> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGlobalIndexer ( ) const
inline

Definition at line 410 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const panzer::GlobalIndexer> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getRangeGlobalIndexer ( ) const
inlinevirtual

Get the range global indexer object associated with this factory.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 413 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const panzer::GlobalIndexer> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getDomainGlobalIndexer ( ) const
inlinevirtual

Get the domain global indexer object associated with this factory.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 416 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
const std::vector<Teuchos::RCP<const GlobalIndexer> >& panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getRangeGlobalIndexers ( ) const
inline

Get global indexers associated with the blocks.

Definition at line 420 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
const std::vector<Teuchos::RCP<const GlobalIndexer> >& panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getDomainGlobalIndexers ( ) const
inline

Get global indexers associated with the blocks.

Definition at line 424 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::addExcludedPair ( int  rowBlock,
int  colBlock 
)

exclude a block pair from the matrix

Definition at line 717 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::addExcludedPairs ( const std::vector< std::pair< int, int > > &  exPairs)

exclude a vector of pairs from the matrix

Definition at line 724 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::initializeContainer_internal ( int  mem,
ThyraObjContainer< double > &  loc 
) const
protected

Initialize container with a specific set of member values.

Note
This will overwrite everything in the container and zero out values not requested.

Definition at line 675 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::initializeGhostedContainer_internal ( int  mem,
ThyraObjContainer< double > &  loc 
) const
protected

Initialize container with a specific set of member values.

Note
This will overwrite everything in the container and zero out values not requested.

Definition at line 696 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const GlobalIndexer > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getGlobalIndexer ( int  i) const
protected
template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const GlobalIndexer > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::getColGlobalIndexer ( int  i) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::makeRoomForBlocks ( std::size_t  blockCnt,
std::size_t  colBlockCnt = 0 
)
protected

Allocate the space in the std::vector objects so we can fill with appropriate Epetra data.

Definition at line 752 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostToGlobalThyraVector ( const Teuchos::RCP< const Thyra::VectorBase< double > > &  in,
const Teuchos::RCP< Thyra::VectorBase< double > > &  out,
bool  col 
) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostToGlobalThyraMatrix ( const Thyra::LinearOpBase< double > &  in,
Thyra::LinearOpBase< double > &  out 
) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::globalToGhostThyraVector ( const Teuchos::RCP< const Thyra::VectorBase< double > > &  in,
const Teuchos::RCP< Thyra::VectorBase< double > > &  out,
bool  col 
) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::adjustForDirichletConditions ( const Epetra_Vector local_bcs,
const Epetra_Vector global_bcs,
const Teuchos::Ptr< Epetra_Vector > &  f,
const Teuchos::Ptr< Epetra_CrsMatrix > &  A,
bool  zeroVectorRows 
) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostToGlobalEpetraVector ( int  i,
const Epetra_Vector in,
Epetra_Vector out,
bool  col 
) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::globalToGhostEpetraVector ( int  i,
const Epetra_Vector in,
Epetra_Vector out,
bool  col 
) const
protected
template<typename Traits , typename LocalOrdinalT >
void panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostToGlobalEpetraMatrix ( int  blockRow,
const Epetra_CrsMatrix in,
Epetra_CrsMatrix out 
) const
protected
template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildMap ( int  i) const
protectedvirtual

Build the i-th owned map from the owned indices of the i-th global indexer.

Parameters
[in]iThe index into the list of global indexers.
Returns
The i-th owned map.

Definition at line 1476 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGhostedMap ( int  i) const
protectedvirtual

Build the i-th ghosted map from the owned and ghosted indices of the i-th global indexer.

Parameters
[in]iThe index into the list of global indexers.
Returns
The i-th owned and ghosted map.

Definition at line 1509 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGhostedMap2 ( int  i) const
protectedvirtual

Build the i-th ghosted map from the ghosted indices of the i-th global indexer.

Parameters
[in]iThe index into the list of global indexers.
Returns
The i-th ghosted map.

Definition at line 1527 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildColMap ( int  i) const
protectedvirtual

Build the i-th owned column map from the owned indices of the i-th (column) global indexer.

Parameters
[in]iThe index into the list of (column) global indexers.
Returns
The i-th owned column map.

Definition at line 1488 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildColGhostedMap ( int  i) const
protectedvirtual

Build the i-th ghosted column map from the owned and ghosted indices of the i-th (column) global indexer.

Parameters
[in]iThe index into the list of (column) global indexers.
Returns
The i-th owned and ghosted column map.

Definition at line 1545 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_Map > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildColGhostedMap2 ( int  i) const
protectedvirtual

Build the i-th ghosted column map from the ghosted indices of the i-th (column) global indexer.

Parameters
[in]iThe index into the list of (column) global indexers.
Returns
The i-th ghosted column map.

Definition at line 1565 of file Panzer_BlockedEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_CrsGraph > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGraph ( int  i,
int  j 
) const
protectedvirtual
template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_CrsGraph > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildGhostedGraph ( int  i,
int  j,
bool  optimizeStorage 
) const
protectedvirtual
template<typename Traits , typename LocalOrdinalT >
const Teuchos::RCP< Epetra_CrsGraph > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::buildFilteredGhostedGraph ( int  i,
int  j 
) const
protectedvirtual

Member Data Documentation

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const DOFManagerContainer> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::rowDOFManagerContainer_
protected

Definition at line 536 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const DOFManagerContainer> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colDOFManagerContainer_
protected

Definition at line 537 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
bool panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::useColGidProviders_
protected

Definition at line 539 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::unordered_set<std::pair<int,int>,panzer::pair_hash> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::excludedPairs_
protected

Definition at line 542 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const Thyra::VectorSpaceBase<double> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::rangeSpace_
mutableprotected

Definition at line 552 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const Thyra::VectorSpaceBase<double> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::domainSpace_
mutableprotected

Definition at line 553 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const Thyra::VectorSpaceBase<double> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostedRangeSpace_
mutableprotected

Definition at line 555 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const Thyra::VectorSpaceBase<double> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostedDomainSpace_
mutableprotected

Definition at line 556 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const Epetra_Comm> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::eComm_
protected

Definition at line 652 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::rawMpiComm_
protected

Definition at line 653 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
Teuchos::RCP<Teuchos::MpiComm<int> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::tComm_
protected

Definition at line 654 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Map> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::maps_
mutableprotected

The list of owned maps corresponding to the owned indices of the global indexers.

Definition at line 660 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Map> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostedMaps_
mutableprotected

The list of ghosted maps corresponding to the owned and ghosted indices of the global indexers.

Definition at line 666 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Map> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostedMaps2_
mutableprotected

The list of ghosted maps corresponding to the ghosted indices of the global indexers.

Definition at line 672 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Import> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::importers_
mutableprotected

The list of ghosted importers corresponding to ghostedMaps_.

Definition at line 677 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Import> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::importers2_
mutableprotected

The list of ghosted importers corresponding to ghostedMaps2_.

Definition at line 682 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Export> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::exporters_
mutableprotected

Definition at line 684 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Map> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colMaps_
mutableprotected

The list of owned column maps corresponding to the owned indices of the (column) global indexers.

Definition at line 690 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Map> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colGhostedMaps_
mutableprotected

The list of ghosted column maps corresponding to the owned and ghosted indices of the (column) global indexers.

Definition at line 696 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Map> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colGhostedMaps2_
mutableprotected

The list of ghosted column maps corresponding to the ghosted indices of the (column) global indexers.

Definition at line 702 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Import> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colImporters_
mutableprotected

The list of ghosted importers corresponding to colGhostedMaps_.

Definition at line 707 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Import> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colImporters2_
mutableprotected

The list of ghosted importers corresponding to colGhostedMaps2_.

Definition at line 713 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::vector<Teuchos::RCP<Epetra_Export> > panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::colExporters_
mutableprotected

Definition at line 715 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::unordered_map<std::pair<int,int>,Teuchos::RCP<Epetra_CrsGraph>,panzer::pair_hash> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::graphs_
mutableprotected

Definition at line 717 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
std::unordered_map<std::pair<int,int>,Teuchos::RCP<Epetra_CrsGraph>,panzer::pair_hash> panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostedGraphs_
mutableprotected

Definition at line 718 of file Panzer_BlockedEpetraLinearObjFactory.hpp.

template<typename Traits, typename LocalOrdinalT>
bool panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >::useDiscreteAdjoint_
protected

Definition at line 720 of file Panzer_BlockedEpetraLinearObjFactory.hpp.


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