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

#include <Panzer_SGEpetraLinearObjFactory_decl.hpp>

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

Public Member Functions

 SGEpetraLinearObjFactory (const Teuchos::RCP< EpetraLinearObjFactory< Traits, LocalOrdinalT > > &epetraFact, const Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > &expansion, const Teuchos::RCP< const EpetraExt::MultiComm > &globalMultiComm)
 
virtual ~SGEpetraLinearObjFactory ()
 
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 mem) const
 
virtual void ghostToGlobalContainer (const LinearObjContainer &ghostContainer, LinearObjContainer &container, int mem) 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::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
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 mem, LinearObjContainer &loc) const
 
void initializeGhostedContainer (int mem, LinearObjContainer &loc) const
 
Teuchos::RCP
< EpetraLinearObjFactory
< Traits, LocalOrdinalT > > 
getEpetraFactory () const
 
Teuchos::RCP< const
Stokhos::OrthogPolyExpansion
< int, double > > 
getExpansion () const
 Accessor for the expansion object. More...
 
Teuchos::RCP
< Stokhos::OrthogPolyExpansion
< int, double > > 
getExpansion ()
 Accessor for the expansion object. More...
 
Teuchos::RCP
< Stokhos::EpetraVectorOrthogPoly > 
getVectorOrthogPoly () const
 Set orthog poly object, this serves as a template for converting vectors to block vectors. More...
 
Teuchos::RCP< const Epetra_MapgetMap ()
 get the map from the matrix, this is the map for the solution vector More...
 
Teuchos::RCP< const Epetra_MapgetSGBlockMap () const
 get the block map needed by Stokhos to describe the parallel layout of the SG unknowns More...
 
Teuchos::RCP< const
panzer::UniqueGlobalIndexerBase
getUniqueGlobalIndexerBase () const
 Get the unique global indexer this factory was created with. 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 &loc) const
 
virtual void endFill (LinearObjContainer &loc) const
 

Protected Attributes

Teuchos::RCP
< EpetraLinearObjFactory
< Traits, LocalOrdinalT > > 
epetraFact_
 
Teuchos::RCP
< Stokhos::OrthogPolyExpansion
< int, double > > 
expansion_
 
Teuchos::RCP< const
EpetraExt::MultiComm
globalMultiComm_
 
Teuchos::RCP< const Epetra_MapsgBlockMap_
 

Detailed Description

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

Linear object factory for constructing Stochastic Galerkin epetra linear object containers. Also handles some of the global to ghosting (and vice-versa) communication.

Definition at line 78 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

Constructor & Destructor Documentation

template<typename Traits , typename LocalOrdinalT >
panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::SGEpetraLinearObjFactory ( const Teuchos::RCP< EpetraLinearObjFactory< Traits, LocalOrdinalT > > &  epetraFact,
const Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > &  expansion,
const Teuchos::RCP< const EpetraExt::MultiComm > &  globalMultiComm 
)

Definition at line 55 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::~SGEpetraLinearObjFactory ( )
virtual

Definition at line 67 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

Member Function Documentation

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

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 89 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

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

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 92 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< LinearObjContainer > panzer::SGEpetraLinearObjFactory< 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 74 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

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

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 89 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< LinearObjContainer > panzer::SGEpetraLinearObjFactory< 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 97 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

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

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 112 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::globalToGhostContainer ( const LinearObjContainer container,
LinearObjContainer ghostContainer,
int  mem 
) const
virtual
template<typename Traits , typename LocalOrdinalT >
void panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::ghostToGlobalContainer ( const LinearObjContainer ghostContainer,
LinearObjContainer container,
int  mem 
) const
virtual
template<typename Traits , typename LocalOrdinalT >
void panzer::SGEpetraLinearObjFactory< 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 179 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::SGEpetraLinearObjFactory< 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 207 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

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

Acess to the MPI Comm used in constructing this LOF.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 215 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

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

Use preconstructed scatter evaluators.

Definition at line 129 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

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

Use preconstructed gather evaluators.

Definition at line 134 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

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

Use preconstructed gather evaluators.

Definition at line 139 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

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

Use preconstructed gather evaluators.

Definition at line 144 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

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

Use preconstructed dirichlet scatter evaluators.

Definition at line 149 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::SGEpetraLinearObjFactory< 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 222 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
void panzer::SGEpetraLinearObjFactory< 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 247 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<EpetraLinearObjFactory<Traits,LocalOrdinalT> > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getEpetraFactory ( ) const
inline

Extract underlying epetra factory from SGEpetraLinearOpFactory

Definition at line 170 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<const Stokhos::OrthogPolyExpansion<int,double> > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getExpansion ( ) const
inline

Accessor for the expansion object.

Definition at line 174 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getExpansion ( )
inline

Accessor for the expansion object.

Definition at line 178 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getVectorOrthogPoly ( ) const

Set orthog poly object, this serves as a template for converting vectors to block vectors.

Definition at line 271 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Epetra_Map > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getMap ( )

get the map from the matrix, this is the map for the solution vector

Definition at line 280 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP< const Epetra_Map > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getSGBlockMap ( ) const

get the block map needed by Stokhos to describe the parallel layout of the SG unknowns

Definition at line 289 of file Panzer_SGEpetraLinearObjFactory_impl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<const panzer::UniqueGlobalIndexerBase> panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::getUniqueGlobalIndexerBase ( ) const
inlinevirtual

Get the unique global indexer this factory was created with.

Implements panzer::LinearObjFactory< Traits >.

Definition at line 191 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

Member Data Documentation

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<EpetraLinearObjFactory<Traits,LocalOrdinalT> > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::epetraFact_
protected

Definition at line 196 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> > panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::expansion_
protected

Definition at line 197 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<const EpetraExt::MultiComm> panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::globalMultiComm_
protected

Definition at line 198 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.

template<typename Traits , typename LocalOrdinalT >
Teuchos::RCP<const Epetra_Map> panzer::SGEpetraLinearObjFactory< Traits, LocalOrdinalT >::sgBlockMap_
mutableprotected

Definition at line 200 of file Panzer_SGEpetraLinearObjFactory_decl.hpp.


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