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

An Epetra operator representing the block stochastic Galerkin operator generated by fully assembling the matrix. The ordering of this operator is interlaced. That means that all stochastic degrees of freedom associated with a deterministic degree of freedom are interlaced. The result is a large sparse matrix that is composed of small (relatively) dense blocks. More...

#include <Stokhos_InterlacedOperator.hpp>

Inheritance diagram for Stokhos::InterlacedOperator:
Inheritance graph
[legend]

Public Member Functions

 InterlacedOperator (const Teuchos::RCP< const EpetraExt::MultiComm > &sg_comm, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &sg_basis, const Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > &epetraCijk, const Teuchos::RCP< const Epetra_CrsGraph > &base_graph, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Constructor. More...
 
virtual ~InterlacedOperator ()
 Destructor. More...
 
void SumIntoGlobalBlock_Deterministic (double alpha, const Epetra_RowMatrix &determBlock, int Row, int Col)
 Sum into global matrix. More...
 
- Public Member Functions inherited from Stokhos::SGOperator
 SGOperator ()
 Constructor. More...
 
virtual ~SGOperator ()
 Destructor. More...
 
- Public Member Functions inherited from Epetra_Operator
virtual int SetUseTranspose (bool UseTranspose)=0
 
virtual int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
 
virtual int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
 
virtual double NormInf () const =0
 
virtual const char * Label () const =0
 
virtual bool UseTranspose () const =0
 
virtual bool HasNormInf () const =0
 
virtual const Epetra_CommComm () const =0
 
virtual const Epetra_MapOperatorDomainMap () const =0
 
virtual const Epetra_MapOperatorRangeMap () const =0
 

Protected Types

typedef Stokhos::Sparse3Tensor
< int, double
Cijk_type
 Short-hand for Cijk. More...
 

Protected Attributes

Teuchos::RCP< const
EpetraExt::MultiComm > 
sg_comm
 Stores SG parallel communicator. More...
 
Teuchos::RCP< const
Stokhos::OrthogPolyBasis< int,
double > > 
sg_basis
 Stochastic Galerking basis. More...
 
Teuchos::RCP< const
Stokhos::EpetraSparse3Tensor
epetraCijk
 Stores Epetra Cijk tensor. More...
 
Teuchos::RCP< const Cijk_typeCijk
 Stores triple product tensor. More...
 
Teuchos::RCP
< Stokhos::EpetraOperatorOrthogPoly
block_ops
 Stores operators. More...
 
bool scale_op
 Flag indicating whether operator be scaled with <^2> More...
 
bool include_mean
 Flag indicating whether to include mean term. More...
 
bool only_use_linear
 Flag indicating whether to only use linear terms. More...
 
int determOffset_
 

Private Member Functions

 InterlacedOperator (const InterlacedOperator &)
 Private to prohibit copying. More...
 
InterlacedOperatoroperator= (const InterlacedOperator &)
 Private to prohibit copying. More...
 

Stokhos::SGOperator methods

virtual void setupOperator (const Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > &poly)
 Setup operator. More...
 
virtual Teuchos::RCP
< Stokhos::EpetraOperatorOrthogPoly
getSGPolynomial ()
 Get SG polynomial. More...
 
virtual Teuchos::RCP< const
Stokhos::EpetraOperatorOrthogPoly
getSGPolynomial () const
 Get SG polynomial. More...
 

Detailed Description

An Epetra operator representing the block stochastic Galerkin operator generated by fully assembling the matrix. The ordering of this operator is interlaced. That means that all stochastic degrees of freedom associated with a deterministic degree of freedom are interlaced. The result is a large sparse matrix that is composed of small (relatively) dense blocks.

Definition at line 30 of file Stokhos_InterlacedOperator.hpp.

Member Typedef Documentation

Short-hand for Cijk.

Definition at line 87 of file Stokhos_InterlacedOperator.hpp.

Constructor & Destructor Documentation

Stokhos::InterlacedOperator::InterlacedOperator ( const Teuchos::RCP< const EpetraExt::MultiComm > &  sg_comm,
const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &  sg_basis,
const Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > &  epetraCijk,
const Teuchos::RCP< const Epetra_CrsGraph > &  base_graph,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)

Constructor.

Definition at line 17 of file Stokhos_InterlacedOperator.cpp.

Stokhos::InterlacedOperator::~InterlacedOperator ( )
virtual

Destructor.

Definition at line 42 of file Stokhos_InterlacedOperator.cpp.

Stokhos::InterlacedOperator::InterlacedOperator ( const InterlacedOperator )
private

Private to prohibit copying.

Member Function Documentation

void Stokhos::InterlacedOperator::setupOperator ( const Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > &  poly)
virtual

Setup operator.

Implements Stokhos::SGOperator.

Definition at line 48 of file Stokhos_InterlacedOperator.cpp.

Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > Stokhos::InterlacedOperator::getSGPolynomial ( )
virtual

Get SG polynomial.

Implements Stokhos::SGOperator.

Definition at line 90 of file Stokhos_InterlacedOperator.cpp.

Teuchos::RCP< const Stokhos::EpetraOperatorOrthogPoly > Stokhos::InterlacedOperator::getSGPolynomial ( ) const
virtual

Get SG polynomial.

Implements Stokhos::SGOperator.

Definition at line 97 of file Stokhos_InterlacedOperator.cpp.

void Stokhos::InterlacedOperator::SumIntoGlobalBlock_Deterministic ( double  alpha,
const Epetra_RowMatrix determBlock,
int  Row,
int  Col 
)

Sum into global matrix.

Definition at line 103 of file Stokhos_InterlacedOperator.cpp.

InterlacedOperator& Stokhos::InterlacedOperator::operator= ( const InterlacedOperator )
private

Private to prohibit copying.

Member Data Documentation

Teuchos::RCP<const EpetraExt::MultiComm> Stokhos::InterlacedOperator::sg_comm
protected

Stores SG parallel communicator.

Definition at line 78 of file Stokhos_InterlacedOperator.hpp.

Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> > Stokhos::InterlacedOperator::sg_basis
protected

Stochastic Galerking basis.

Definition at line 81 of file Stokhos_InterlacedOperator.hpp.

Teuchos::RCP<const Stokhos::EpetraSparse3Tensor> Stokhos::InterlacedOperator::epetraCijk
protected

Stores Epetra Cijk tensor.

Definition at line 84 of file Stokhos_InterlacedOperator.hpp.

Teuchos::RCP<const Cijk_type> Stokhos::InterlacedOperator::Cijk
protected

Stores triple product tensor.

Definition at line 90 of file Stokhos_InterlacedOperator.hpp.

Teuchos::RCP<Stokhos::EpetraOperatorOrthogPoly > Stokhos::InterlacedOperator::block_ops
protected

Stores operators.

Definition at line 93 of file Stokhos_InterlacedOperator.hpp.

bool Stokhos::InterlacedOperator::scale_op
protected

Flag indicating whether operator be scaled with <^2>

Definition at line 96 of file Stokhos_InterlacedOperator.hpp.

bool Stokhos::InterlacedOperator::include_mean
protected

Flag indicating whether to include mean term.

Definition at line 99 of file Stokhos_InterlacedOperator.hpp.

bool Stokhos::InterlacedOperator::only_use_linear
protected

Flag indicating whether to only use linear terms.

Definition at line 102 of file Stokhos_InterlacedOperator.hpp.

int Stokhos::InterlacedOperator::determOffset_
protected

Definition at line 104 of file Stokhos_InterlacedOperator.hpp.


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