Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Private Types | Private Attributes | List of all members
panzer::LinearObjFactory< Traits > Class Template Referenceabstract

#include <Panzer_LinearObjFactory.hpp>

Inheritance diagram for panzer::LinearObjFactory< Traits >:
Inheritance graph
[legend]

Classes

struct  Gather_Builder
 
struct  GatherDomain_Builder
 
struct  GatherOrientation_Builder
 
struct  GatherTangent_Builder
 
struct  Scatter_Builder
 
struct  ScatterDirichlet_Builder
 

Public Member Functions

virtual ~LinearObjFactory ()
 
template<typename BuilderT >
void buildGatherScatterEvaluators (const BuilderT &builder)
 
virtual void readVector (const std::string &identifier, LinearObjContainer &loc, int id) const =0
 
virtual void writeVector (const std::string &identifier, const LinearObjContainer &loc, int id) const =0
 
virtual Teuchos::RCP
< LinearObjContainer
buildLinearObjContainer () const =0
 
virtual Teuchos::RCP
< LinearObjContainer
buildPrimitiveLinearObjContainer () const =0
 
virtual Teuchos::RCP
< LinearObjContainer
buildGhostedLinearObjContainer () const =0
 
virtual Teuchos::RCP
< LinearObjContainer
buildPrimitiveGhostedLinearObjContainer () const =0
 
virtual Teuchos::RCP
< ReadOnlyVector_GlobalEvaluationData
buildReadOnlyDomainContainer () const =0
 
virtual void globalToGhostContainer (const LinearObjContainer &container, LinearObjContainer &ghostContainer, int) const =0
 
virtual void ghostToGlobalContainer (const LinearObjContainer &ghostContainer, LinearObjContainer &container, int) const =0
 
virtual void initializeContainer (int, LinearObjContainer &loc) const =0
 
virtual void initializeGhostedContainer (int, LinearObjContainer &loc) const =0
 
virtual void adjustForDirichletConditions (const LinearObjContainer &localBCRows, const LinearObjContainer &globalBCRows, LinearObjContainer &ghostedObjs, bool zeroVectorRows=false, bool adjustX=false) const =0
 
virtual void applyDirichletBCs (const LinearObjContainer &counter, LinearObjContainer &result) const =0
 
virtual Teuchos::MpiComm< int > getComm () const =0
 
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 Teuchos::RCP< const
panzer::GlobalIndexer
getDomainGlobalIndexer () const =0
 Get the domain global indexer object associated with this factory. More...
 
virtual Teuchos::RCP< const
panzer::GlobalIndexer
getRangeGlobalIndexer () const =0
 Get the range global indexer object associated with this factory. More...
 
virtual void beginFill (LinearObjContainer &) const
 
virtual void endFill (LinearObjContainer &) const
 

Private Types

typedef PHX::TemplateManager
< typename Traits::EvalTypes,
panzer::CloneableEvaluator,
PHX::EvaluatorDerived< _,
Traits > > 
Evaluator_TemplateManager
 

Private Attributes

Teuchos::RCP
< Evaluator_TemplateManager
scatterManager_
 
Teuchos::RCP
< Evaluator_TemplateManager
scatterDirichletManager_
 
Teuchos::RCP
< Evaluator_TemplateManager
gatherManager_
 
Teuchos::RCP
< Evaluator_TemplateManager
gatherTangentManager_
 
Teuchos::RCP
< Evaluator_TemplateManager
gatherDomainManager_
 
Teuchos::RCP
< Evaluator_TemplateManager
gatherOrientManager_
 

Detailed Description

template<typename Traits>
class panzer::LinearObjFactory< Traits >

Abstract factory that builds the linear algebra objects required for the assembly including the gather/scatter evaluator objects.

The interface for construction of the gather scatter is externally very simple, but in fact under the hood it is quite complex. The user of this factory object simply calls buildGather(const Teuchos::ParameterList & pl) const, buildScatter(const Teuchos::ParameterList & pl) const, or buildScatterDirichlet(const Teuchos::ParameterList & pl) const.

To implement a version of this class an author must overide all the linear algebra construction functions. The new version should also call the base class version of buildGatherScatterEvaluators in the constructor with a reference to itself passed in as an argument. This requires the new version of the class to implement the following functions

This builds the correct scatter/gather/scatter-dirichlet evaluator objects and returns them as a CloneableEvaluator (These evaluators must overide the CloneableEvaluator::clone function which takes a parameter list). The cloned evaluators will be the ones actually returned from the buildGather(const Teuchos::ParameterList & pl) const, buildGatherDomain(const Teuchos::ParameterList & pl) const, buildGatherOrientation(const Teuchos::ParameterList & pl) const, buildScatter(const Teuchos::ParameterList & pl) const, or buildScatterDirichlet(const Teuchos::ParameterList & pl) const functions.

Definition at line 104 of file Panzer_LinearObjFactory.hpp.

Member Typedef Documentation

Definition at line 292 of file Panzer_LinearObjFactory.hpp.

Constructor & Destructor Documentation

template<typename Traits>
virtual panzer::LinearObjFactory< Traits >::~LinearObjFactory ( )
inlinevirtual

Definition at line 106 of file Panzer_LinearObjFactory.hpp.

Member Function Documentation

template<typename Traits >
template<typename BuilderT >
void panzer::LinearObjFactory< Traits >::buildGatherScatterEvaluators ( const BuilderT &  builder)
inline

This builds all the required evaluators. It is required to be called before the build[Gather,Scatter,ScatterDirichlet] functions are called. This would typically be called by the inheriting class.

Parameters
[in]builderTemplate class to build all required evaluators. The class has the following interface.

Definition at line 372 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::readVector ( const std::string &  identifier,
LinearObjContainer loc,
int  id 
) const
pure 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

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::writeVector ( const std::string &  identifier,
const LinearObjContainer loc,
int  id 
) const
pure 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

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual Teuchos::RCP<LinearObjContainer> panzer::LinearObjFactory< Traits >::buildLinearObjContainer ( ) const
pure virtual
template<typename Traits>
virtual Teuchos::RCP<LinearObjContainer> panzer::LinearObjFactory< Traits >::buildPrimitiveLinearObjContainer ( ) const
pure 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.

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual Teuchos::RCP<LinearObjContainer> panzer::LinearObjFactory< Traits >::buildGhostedLinearObjContainer ( ) const
pure virtual
template<typename Traits>
virtual Teuchos::RCP<LinearObjContainer> panzer::LinearObjFactory< Traits >::buildPrimitiveGhostedLinearObjContainer ( ) const
pure 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.

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual Teuchos::RCP<ReadOnlyVector_GlobalEvaluationData> panzer::LinearObjFactory< Traits >::buildReadOnlyDomainContainer ( ) const
pure 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.

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::globalToGhostContainer ( const LinearObjContainer container,
LinearObjContainer ghostContainer,
int   
) const
pure virtual
template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::ghostToGlobalContainer ( const LinearObjContainer ghostContainer,
LinearObjContainer container,
int   
) const
pure virtual
template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::initializeContainer ( int  ,
LinearObjContainer loc 
) const
pure virtual

Initialize container with a specific set of member values.

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

Implemented in panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, and panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::initializeGhostedContainer ( int  ,
LinearObjContainer loc 
) const
pure virtual

Initialize container with a specific set of member values.

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

Implemented in panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, and panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::adjustForDirichletConditions ( const LinearObjContainer localBCRows,
const LinearObjContainer globalBCRows,
LinearObjContainer ghostedObjs,
bool  zeroVectorRows = false,
bool  adjustX = false 
) const
pure 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.

Parameters
[in]localBCRowsLinear object container uses the X vector to indicate locally set dirichlet conditions. The format is if an entry of the vector is nonzero then it was set as a dirichlet condition.
[in]globalBCRowsLinear object container uses the X vector to indicate globally set dirichlet conditions. The format is if an entry of the vector is nonzero then it was set as a dirichlet condition.
[in,out]ghostedObjsGhosted linear object container storing the residual and jacobian matrix for any boundary conditions set. The matrix will be modified by zeroing any rows that are set as nonzero in globalBCRows but zero in localBCRows (similarly for the vector). If a row is nonzero in both localBCRows and globalBCRows then those rows in both the matrix and the residual vector are devided by the corresponding entry in the globalBCRows.
[in]zeroVectorRowsInstead of preserving (and scaling) the vector rows, setting this to true will zero them instead.

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

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

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

Parameters
[in]counterContains a counter vector (the "x" vector) indicating which rows contain dirichlet conditions by having a non zero entry. The dirichlet condition values are specified in the "f" vector.
[in,out]resultThe vector to be modifed is the "f" vector.

Implemented in panzer::BlockedTpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >, panzer::BlockedEpetraLinearObjFactory< Traits, LocalOrdinalT >, and panzer::TpetraLinearObjFactory< Traits, ScalarT, LocalOrdinalT, GlobalOrdinalT, NodeT >.

template<typename Traits>
virtual Teuchos::MpiComm<int> panzer::LinearObjFactory< Traits >::getComm ( ) const
pure virtual
template<typename Traits>
template<typename EvalT >
Teuchos::RCP<PHX::Evaluator<Traits> > panzer::LinearObjFactory< Traits >::buildScatter ( const Teuchos::ParameterList pl) const
inline

Use preconstructed scatter evaluators.

Definition at line 251 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
template<typename EvalT >
Teuchos::RCP<PHX::Evaluator<Traits> > panzer::LinearObjFactory< Traits >::buildGather ( const Teuchos::ParameterList pl) const
inline

Use preconstructed gather evaluators.

Definition at line 256 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
template<typename EvalT >
Teuchos::RCP<PHX::Evaluator<Traits> > panzer::LinearObjFactory< Traits >::buildGatherTangent ( const Teuchos::ParameterList pl) const
inline

Use preconstructed gather evaluators.

Definition at line 261 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
template<typename EvalT >
Teuchos::RCP<PHX::Evaluator<Traits> > panzer::LinearObjFactory< Traits >::buildGatherDomain ( const Teuchos::ParameterList pl) const
inline

Use preconstructed gather evaluators.

Definition at line 266 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
template<typename EvalT >
Teuchos::RCP<PHX::Evaluator<Traits> > panzer::LinearObjFactory< Traits >::buildGatherOrientation ( const Teuchos::ParameterList pl) const
inline

Use preconstructed gather evaluators.

Definition at line 271 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
template<typename EvalT >
Teuchos::RCP<PHX::Evaluator<Traits> > panzer::LinearObjFactory< Traits >::buildScatterDirichlet ( const Teuchos::ParameterList pl) const
inline

Use preconstructed dirichlet scatter evaluators.

Definition at line 276 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
virtual Teuchos::RCP<const panzer::GlobalIndexer> panzer::LinearObjFactory< Traits >::getDomainGlobalIndexer ( ) const
pure virtual
template<typename Traits>
virtual Teuchos::RCP<const panzer::GlobalIndexer> panzer::LinearObjFactory< Traits >::getRangeGlobalIndexer ( ) const
pure virtual
template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::beginFill ( LinearObjContainer ) const
inlinevirtual
template<typename Traits>
virtual void panzer::LinearObjFactory< Traits >::endFill ( LinearObjContainer ) const
inlinevirtual

Member Data Documentation

template<typename Traits>
Teuchos::RCP<Evaluator_TemplateManager> panzer::LinearObjFactory< Traits >::scatterManager_
private

Definition at line 295 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
Teuchos::RCP<Evaluator_TemplateManager> panzer::LinearObjFactory< Traits >::scatterDirichletManager_
private

Definition at line 296 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
Teuchos::RCP<Evaluator_TemplateManager> panzer::LinearObjFactory< Traits >::gatherManager_
private

Definition at line 297 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
Teuchos::RCP<Evaluator_TemplateManager> panzer::LinearObjFactory< Traits >::gatherTangentManager_
private

Definition at line 298 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
Teuchos::RCP<Evaluator_TemplateManager> panzer::LinearObjFactory< Traits >::gatherDomainManager_
private

Definition at line 299 of file Panzer_LinearObjFactory.hpp.

template<typename Traits>
Teuchos::RCP<Evaluator_TemplateManager> panzer::LinearObjFactory< Traits >::gatherOrientManager_
private

Definition at line 300 of file Panzer_LinearObjFactory.hpp.


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