LifeV
MonolithicBlock Class Referenceabstract

MonolithicBlock - This is a pure virtual class for the linear operators with a block structure. More...

#include <MonolithicBlock.hpp>

+ Inheritance diagram for MonolithicBlock:
+ Collaboration diagram for MonolithicBlock:

Public Member Functions

virtual void replace_bch (bchandlerPtr_Type &, UInt)
 replaces a BCHandler More...
 
virtual void applyBoundaryConditions (const Real &time)
 Applies the correspondent boundary conditions to every block. More...
 
virtual void applyBoundaryConditions (const Real &time, const UInt block)
 Applies the correspondent boundary conditions to a specified block. More...
 
virtual void resetBlocks ()
 resets the blocks (frees the shared pointers) More...
 
virtual void setComm (std::shared_ptr< Epetra_Comm > comm)
 sets the communicator More...
 
virtual void reset ()
 resets the blocks, boundary conditions, FE spaces. More...
 
virtual void setRobin (matrixPtr_Type &, vectorPtr_Type &)
 Applies the robin preconditioners. More...
 
void couplingMatrix (matrixPtr_Type &bigMatrix, Int flag, const std::vector< fespacePtr_Type > &problem, const std::vector< UInt > &offset, const std::map< ID, ID > &locDofMap, const vectorPtr_Type &numerationInterface, const Real &timeStep=1.e-3, const Real &value=1., const Real &coefficient=1., const Real &rescaleFactor=1.)
 builds the coupling matrix. More...
 
void setConditions (std::vector< bchandlerPtr_Type > &vec)
 sets the vector of raw pointer to the BCHandler More...
 
void setSpaces (std::vector< fespacePtr_Type > &vec)
 sets the vector of raw pointer to the FESpaces More...
 
void setOffsets (UInt blocks,...)
 sets the vector of raw pointer to the offsets of the different blocks More...
 
void robinCoupling (MonolithicBlock::matrixPtr_Type &matrix, Real &alphaf, Real &alphas, UInt coupling, const MonolithicBlock::fespacePtr_Type &FESpace1, const UInt &offset1, const MonolithicBlock::fespacePtr_Type &FESpace2, const UInt &offset2, const std::map< ID, ID > &locDofMap, const MonolithicBlock::vectorPtr_Type &numerationInterface)
 computes the Robin coupling matrix More...
 
virtual void addToBlock (const matrixPtr_Type &Mat, UInt position)
 
virtual void push_back_oper (MonolithicBlock &Oper)
 Pushes a new block. More...
 

Public Types

typedef VectorEpetra vector_Type
 
typedef std::shared_ptr< vector_TypevectorPtr_Type
 
typedef MatrixEpetra< Realmatrix_Type
 
typedef std::shared_ptr< matrix_TypematrixPtr_Type
 
typedef std::shared_ptr< Epetra_OperatorepetraOperatorPtr_Type
 
typedef std::shared_ptr< Preconditionerepetra_preconditioner_ptrtype
 
typedef matrix_Type::matrix_type epetraMatrix_Type
 
typedef SolverAztecOO solver_Type
 
typedef std::shared_ptr< SolverAztecOOsolverPtr_Type
 
typedef FESpace< RegionMesh< LinearTetra >, MapEpetrafespace_Type
 
typedef std::shared_ptr< fespace_TypefespacePtr_Type
 
typedef std::shared_ptr< MapEpetramapPtr_Type
 
typedef std::shared_ptr< BCHandlerbchandlerPtr_Type
 

Constructor & Destructor

 MonolithicBlock ()
 Empty Constructor. More...
 
virtual ~MonolithicBlock ()
 Destructor. More...
 

Pure virtual methods

virtual int solveSystem (const vector_Type &rhs, vector_Type &result, solverPtr_Type &linearSolver)=0
 Solves the preconditioned linear system (used only when dealing with a preconditioner) More...
 
virtual void setDataFromGetPot (const GetPot &data, const std::string &section)=0
 Sets the parameters needed by the preconditioner from data file. More...
 
virtual void setupSolver (solver_Type &, const GetPot &)
 Sets the parameters needed by the preconditioner from data file. More...
 
virtual void push_back_matrix (const matrixPtr_Type &Mat, const bool recompute)=0
 pushes a block at the end of the vector More...
 
virtual void addToCoupling (const matrixPtr_Type &Mat, UInt position)=0
 
virtual void addToCoupling (const Real &entry, UInt row, UInt col, UInt position)=0
 
virtual void setRecompute (UInt, bool)
 If not present in the derived class it must not be called (gives an assertion fail) More...
 
virtual void replace_matrix (const matrixPtr_Type &Mat, UInt index)=0
 replaces a block More...
 
virtual void replace_coupling (const matrixPtr_Type &Mat, UInt index)=0
 replaces a coupling block More...
 
virtual void GlobalAssemble ()=0
 runs GlobalAssemble on the blocks More...
 
virtual void blockAssembling ()
 sums the coupling matrices with the corresponding blocks More...
 
virtual void coupler (mapPtr_Type &map, const std::map< ID, ID > &locDofMap, const vectorPtr_Type &numerationInterface, const Real &timeStep, const Real &coefficient, const Real &rescaleFactor)=0
 Computes the coupling. More...
 
virtual void coupler (mapPtr_Type &map, const std::map< ID, ID > &locDofMap, const vectorPtr_Type &numerationInterface, const Real &timeStep, const Real &coefficient, const Real &rescaleFactor, UInt couplingFlag)=0
 Adds a default coupling block at a specified position. More...
 
virtual bool set ()=0
 returns true if the operator is set More...
 
virtual void push_back_coupling (matrixPtr_Type &coupling)=0
 Pushes a new coupling matrix in the corresponding vector. More...
 
virtual void replace_precs (const epetraOperatorPtr_Type &, UInt)
 replace a block preconditioner More...
 
virtual void push_back_precs (const epetraOperatorPtr_Type &)
 pushes back a block preconditioner More...
 

Get Methods

const std::vector< matrixPtr_Type > & blockVector ()
 returns the vector of pointers to the blocks (by const reference). More...
 
const std::vector< bchandlerPtr_Type > & BChVector ()
 returns the vector of pointers to the BCHandlers (by const reference). More...
 
const std::vector< fespacePtr_Type > & FESpaceVector ()
 returns the vector of pointers to the FE spaces (by const reference). More...
 
const std::vector< UInt > & offsetVector ()
 returns the vector of the offsets (by const reference). More...
 
virtual const UInt whereIsBlock (UInt position) const =0
 

Protected Methods

virtual void blockAssembling (const UInt)
 sums the coupling matrix in the specified position with the corresponding block More...
 
template<typename Operator >
void swap (std::shared_ptr< Operator > &operFrom, std::shared_ptr< Operator > &OperTo)
 swaps two std::shared_ptr. The tamplate argument of the shared_ptr is templated More...
 
template<typename Operator >
void insert (std::vector< Operator > &operFrom, std::vector< Operator > &OperTo)
 swaps two std::shared_ptr. The tamplate argument of the shared_ptr is templated More...
 

Protected Members

std::vector< bchandlerPtr_TypeM_bch
 
std::vector< matrixPtr_TypeM_blocks
 
std::vector< fespacePtr_TypeM_FESpace
 
std::vector< UIntM_offset
 
vectorPtr_Type M_numerationInterface
 
std::shared_ptr< Epetra_Comm > M_comm
 

Private Methods

 MonolithicBlock (const MonolithicBlock &)
 private copy constructor: this class should not be copied. More...
 

Detailed Description

MonolithicBlock - This is a pure virtual class for the linear operators with a block structure.

(i.e. block matrices and preconditioners).

Author
Paolo Crosetto

The specializations of this class can be used to handle many types of composed preconditioners and composed operators. It conteins a vector of pointers to EpetraMatrix holding the diferent blocks. These matrices should be zero everywhere but on the diagonal block corresponding to the corresponding problem. This class also handles the coupling through a vector of pointer to coupling matrices.

Definition at line 72 of file MonolithicBlock.hpp.

Member Typedef Documentation

◆ vector_Type

Definition at line 78 of file MonolithicBlock.hpp.

◆ vectorPtr_Type

typedef std::shared_ptr< vector_Type > vectorPtr_Type

Definition at line 79 of file MonolithicBlock.hpp.

◆ matrix_Type

Definition at line 80 of file MonolithicBlock.hpp.

◆ matrixPtr_Type

typedef std::shared_ptr< matrix_Type > matrixPtr_Type

Definition at line 81 of file MonolithicBlock.hpp.

◆ epetraOperatorPtr_Type

typedef std::shared_ptr< Epetra_Operator > epetraOperatorPtr_Type

Definition at line 82 of file MonolithicBlock.hpp.

◆ epetra_preconditioner_ptrtype

typedef std::shared_ptr< Preconditioner > epetra_preconditioner_ptrtype

Definition at line 83 of file MonolithicBlock.hpp.

◆ epetraMatrix_Type

◆ solver_Type

Definition at line 85 of file MonolithicBlock.hpp.

◆ solverPtr_Type

typedef std::shared_ptr< SolverAztecOO > solverPtr_Type

Definition at line 86 of file MonolithicBlock.hpp.

◆ fespace_Type

Definition at line 87 of file MonolithicBlock.hpp.

◆ fespacePtr_Type

typedef std::shared_ptr< fespace_Type > fespacePtr_Type

Definition at line 88 of file MonolithicBlock.hpp.

◆ mapPtr_Type

typedef std::shared_ptr< MapEpetra > mapPtr_Type

Definition at line 89 of file MonolithicBlock.hpp.

◆ bchandlerPtr_Type

typedef std::shared_ptr< BCHandler > bchandlerPtr_Type

Definition at line 90 of file MonolithicBlock.hpp.

Constructor & Destructor Documentation

◆ MonolithicBlock() [1/2]

MonolithicBlock ( )
inline

Empty Constructor.

Definition at line 97 of file MonolithicBlock.hpp.

◆ ~MonolithicBlock()

virtual ~MonolithicBlock ( )
inlinevirtual

Destructor.

Definition at line 115 of file MonolithicBlock.hpp.

◆ MonolithicBlock() [2/2]

MonolithicBlock ( const MonolithicBlock )
inlineprivate

private copy constructor: this class should not be copied.

If you need a copy you should implement it, so that it copies the shared pointer one by one, without copying the content.

Definition at line 530 of file MonolithicBlock.hpp.

Member Function Documentation

◆ solveSystem()

virtual int solveSystem ( const vector_Type rhs,
vector_Type result,
solverPtr_Type linearSolver 
)
pure virtual

Solves the preconditioned linear system (used only when dealing with a preconditioner)

Provided the linear solver and the right hand side this method computes the solution and returns it into the result vector.

Parameters
rhsright hand side of the linear system
resultoutput result
linearSolverthe linear system

Implemented in MonolithicBlockMatrix, MonolithicBlockComposedDN, MonolithicBlockComposed, and MonolithicBlockComposedNN.

◆ setDataFromGetPot()

virtual void setDataFromGetPot ( const GetPot data,
const std::string &  section 
)
pure virtual

Sets the parameters needed by the preconditioner from data file.

Parameters
dataGetPot object reading the text data file
sectionstring specifying the path in the data file where to find the options for the operator

Implemented in MonolithicBlockComposed, MonolithicBlockComposedDN, MonolithicBlockMatrix, MonolithicBlockMatrixRN, and MonolithicBlockComposedNN.

◆ setupSolver()

virtual void setupSolver ( solver_Type ,
const GetPot  
)
inlinevirtual

Sets the parameters needed by the preconditioner from data file.

Parameters
dataGetPot object reading the text data file
sectionstring specifying the path in the data file where to find the options for the operator

Reimplemented in MonolithicBlockMatrix.

Definition at line 147 of file MonolithicBlock.hpp.

◆ push_back_matrix()

virtual void push_back_matrix ( const matrixPtr_Type Mat,
const bool  recompute 
)
pure virtual

pushes a block at the end of the vector

adds a new block

Parameters
Matblock matrix to push
recomputeflag stating wether the preconditioner for this block have to be recomputed at every time step

Implemented in MonolithicBlockMatrix, MonolithicBlockComposed, and MonolithicBlockComposedNN.

◆ addToCoupling() [1/2]

virtual void addToCoupling ( const matrixPtr_Type Mat,
UInt  position 
)
pure virtual

adds a new block

Parameters
Matblock matrix to push
recomputeflag stating wether the preconditioner for this block have to be recomputed at every time step

Implemented in MonolithicBlockMatrix, and MonolithicBlockComposed.

◆ addToCoupling() [2/2]

virtual void addToCoupling ( const Real entry,
UInt  row,
UInt  col,
UInt  position 
)
pure virtual

adds an entry to the coupling matrix

Parameters
entryentry
rowrow for the insertion
colcolon for the insertion

Implemented in MonolithicBlockMatrix, and MonolithicBlockComposed.

◆ setRecompute()

virtual void setRecompute ( UInt  ,
bool   
)
inlinevirtual

If not present in the derived class it must not be called (gives an assertion fail)

Reimplemented in MonolithicBlockComposed.

Definition at line 175 of file MonolithicBlock.hpp.

◆ replace_matrix()

virtual void replace_matrix ( const matrixPtr_Type Mat,
UInt  index 
)
pure virtual

replaces a block

replaces a block on a specified position in the vector

Parameters
Matblock matrix to push
indexposition in the vector

Implemented in MonolithicBlockComposed, MonolithicBlockMatrix, and MonolithicBlockComposedNN.

◆ replace_coupling()

virtual void replace_coupling ( const matrixPtr_Type Mat,
UInt  index 
)
pure virtual

replaces a coupling block

replaces a block on a specified position in the vector

Parameters
Matblock matrix to push
indexposition in the vector

Implemented in MonolithicBlockComposed, and MonolithicBlockMatrix.

◆ GlobalAssemble()

virtual void GlobalAssemble ( )
pure virtual

runs GlobalAssemble on the blocks

closes and distributes all the matrices before computing the preconditioner

Implemented in MonolithicBlockComposed, MonolithicBlockMatrix, and MonolithicBlockMatrixRN.

◆ blockAssembling() [1/2]

virtual void blockAssembling ( )
inlinevirtual

sums the coupling matrices with the corresponding blocks

Everything (but the boundary conditions) must have been set before calling this

Reimplemented in MonolithicBlockMatrix, MonolithicBlockComposed, MonolithicBlockMatrixRN, and MonolithicBlockComposedDND.

Definition at line 208 of file MonolithicBlock.hpp.

◆ coupler() [1/2]

virtual void coupler ( mapPtr_Type map,
const std::map< ID, ID > &  locDofMap,
const vectorPtr_Type numerationInterface,
const Real timeStep,
const Real coefficient,
const Real rescaleFactor 
)
pure virtual

Computes the coupling.

computes all the coupling blocks specific for the chosen preconditioner. The coupling is handled through an augmented formulation, introducing new variables (multipliers). Needs as input: the global map of the problem, the FESpaces of the subproblems to be coupled with their offsets, a std::map holding the two numerations of the interface between the two subproblems (the numeration can be different but the nodes must be matching in each subdomain), an EpetraVector defined on the multipliers map containing the corresponding dof at the interface (NB: the multipliers map should be constructed from the second numeration in the std::map). Note that the FESpaces and the offsets have to be set before calling this method.

Parameters
mapthe map of the global problem
locDofMapstd::map with the correspondence between the interface dofs for the two different maps in the subproblems
numerationInterfacevector containing the correspondence of the Lagrange multipliers with the interface dofs

Implemented in MonolithicBlockComposedDN, MonolithicBlockMatrix, MonolithicBlockMatrixRN, and MonolithicBlockComposedNN.

◆ coupler() [2/2]

virtual void coupler ( mapPtr_Type map,
const std::map< ID, ID > &  locDofMap,
const vectorPtr_Type numerationInterface,
const Real timeStep,
const Real coefficient,
const Real rescaleFactor,
UInt  couplingFlag 
)
pure virtual

Adds a default coupling block at a specified position.

The coupling is handled through an augmented formulation, introducing new variables (multipliers). Needs as input: the global map of the problem, the two FESpaces of the subproblems to be coupled with their offsets, a std::map holding the two numerations of the interface between the two subproblems (the numeration can be different but the nodes must be matching in each subdomain), an EpetraVector defined on the multipliers map containing the corresponding dof at the interface (NB: the multipliers map should be constructed from the second numeration in the std::map)

Parameters
mapthe map of the global problem
FESpace1FESpace of the first problem
offset1offset for the first block in the global matrix
FESpace2FESpace of the second problem
offset2offset for the second block in the global matrix
locDofMapstd::map with the correspondence between the interface dofs for the two different maps in the subproblems
numerationInterfacevector containing the correspondence of the Lagrange multipliers with the interface dofs
couplingBlockflag specifying the block associated with the coupling. See MonolithicBlock::couplingMatrix to understand what the values for this flag correspond to.
couplingFlagflag specifying which block must be coupled whith which block.

Implemented in MonolithicBlockComposed, MonolithicBlockMatrix, and MonolithicBlockMatrixRN.

◆ set()

virtual bool set ( )
pure virtual

returns true if the operator is set

Implemented in MonolithicBlockComposedDN, MonolithicBlockMatrix, MonolithicBlockComposedNN, and MonolithicBlockComposed.

◆ push_back_coupling()

virtual void push_back_coupling ( matrixPtr_Type coupling)
pure virtual

Pushes a new coupling matrix in the corresponding vector.

Implemented in MonolithicBlockMatrix, and MonolithicBlockComposed.

◆ replace_precs()

virtual void replace_precs ( const epetraOperatorPtr_Type ,
UInt   
)
inlinevirtual

replace a block preconditioner

(only used if the operator is a preconditioner)

Reimplemented in MonolithicBlockMatrix.

Definition at line 275 of file MonolithicBlock.hpp.

◆ push_back_precs()

virtual void push_back_precs ( const epetraOperatorPtr_Type )
inlinevirtual

pushes back a block preconditioner

(only used if the operator is a preconditioner)

Definition at line 284 of file MonolithicBlock.hpp.

◆ replace_bch()

virtual void replace_bch ( bchandlerPtr_Type ,
UInt   
)
inlinevirtual

replaces a BCHandler

replaces a BCHandler in the vector of bcHandlers at a specified position

Parameters
bchinput BCHandler to replace
positionposition

Definition at line 296 of file MonolithicBlock.hpp.

◆ applyBoundaryConditions() [1/2]

void applyBoundaryConditions ( const Real time)
virtual

Applies the correspondent boundary conditions to every block.

note that this method must be called after blockAssembling(), that sums the coupling conditions to the blocks

Parameters
timetime

Reimplemented in MonolithicBlockMatrix.

Definition at line 117 of file MonolithicBlock.cpp.

◆ applyBoundaryConditions() [2/2]

void applyBoundaryConditions ( const Real time,
const UInt  block 
)
virtual

Applies the correspondent boundary conditions to a specified block.

note that this method must be called after blockAssembling(), that sums the coupling conditions to the blocks

Parameters
timetime
blocknumber of the block

Reimplemented in MonolithicBlockMatrix, and MonolithicBlockComposedNN.

Definition at line 127 of file MonolithicBlock.cpp.

◆ resetBlocks()

virtual void resetBlocks ( )
inlinevirtual

resets the blocks (frees the shared pointers)

Definition at line 317 of file MonolithicBlock.hpp.

◆ setComm()

virtual void setComm ( std::shared_ptr< Epetra_Comm >  comm)
inlinevirtual

sets the communicator

Reimplemented in MonolithicBlockComposedDN.

Definition at line 325 of file MonolithicBlock.hpp.

◆ reset()

virtual void reset ( )
inlinevirtual

resets the blocks, boundary conditions, FE spaces.

Definition at line 333 of file MonolithicBlock.hpp.

◆ setRobin()

virtual void setRobin ( matrixPtr_Type ,
vectorPtr_Type  
)
inlinevirtual

Applies the robin preconditioners.

Applies the robin preconditioners when needed, otherwise does nothing

Reimplemented in MonolithicBlockMatrixRN.

Definition at line 345 of file MonolithicBlock.hpp.

◆ couplingMatrix()

void couplingMatrix ( matrixPtr_Type bigMatrix,
Int  flag,
const std::vector< fespacePtr_Type > &  problem,
const std::vector< UInt > &  offset,
const std::map< ID, ID > &  locDofMap,
const vectorPtr_Type numerationInterface,
const Real timeStep = 1.e-3,
const Real value = 1.,
const Real coefficient = 1.,
const Real rescaleFactor = 1. 
)

builds the coupling matrix.

Public Methods.

Computes the matrix that couples 2 (or more) blocks using an augmented formulation (Lagrange multipliers). It adds from 1 to 4 coupling blocks, depending on the flag 'coupling'. This flag is an Int, it works as the bash command chmod. It ranges from 1 to 15 for two blocks. The values from 15 to 31 account for the presence of a third block: in FSI these three blocks are the fluid block, the structure and the fluid mesh displacement blocks. For 'coupling > 31' the method is called recursively.

Parameters
bigMatrixthe coupling matrix to be built
couplingflag specifying what wether to consider all the coupling or to neglect one part
problem1FESpace of the first block
offset1offset of the first block
problem2FESpace of the second block
offset2offset of the second block
locDofMapstd::map with the correspondence between the interface dofs for the two different maps in the subproblems
numerationInterfacevector containing the correspondence of the Lagrange multipliers with the interface dofs
valuevalue to insert in the coupling blocks

Definition at line 40 of file MonolithicBlock.cpp.

◆ setConditions()

void setConditions ( std::vector< bchandlerPtr_Type > &  vec)

sets the vector of raw pointer to the BCHandler

each entry of the vector correspond to a block in the same position in the vector of matrix pointers M_blocks. A vector of any size can be passed.

Parameters
vecthe vector of BCHandlers

Definition at line 139 of file MonolithicBlock.cpp.

◆ setSpaces()

void setSpaces ( std::vector< fespacePtr_Type > &  vec)

sets the vector of raw pointer to the FESpaces

A vector of arbitrary size can be passed

Parameters
vecthe vector of FESpaces

Definition at line 145 of file MonolithicBlock.cpp.

◆ setOffsets()

void setOffsets ( UInt  blocks,
  ... 
)

sets the vector of raw pointer to the offsets of the different blocks

An arbitrary number of raw pointers can be passed.

Parameters
blockstotal number of input offsets

Definition at line 151 of file MonolithicBlock.cpp.

◆ robinCoupling()

void robinCoupling ( MonolithicBlock::matrixPtr_Type matrix,
Real alphaf,
Real alphas,
UInt  coupling,
const MonolithicBlock::fespacePtr_Type FESpace1,
const UInt offset1,
const MonolithicBlock::fespacePtr_Type FESpace2,
const UInt offset2,
const std::map< ID, ID > &  locDofMap,
const MonolithicBlock::vectorPtr_Type numerationInterface 
)

computes the Robin coupling matrix

Computes a matrix that mixes the coupling conditions between the blocks: [0,0,0,0;0,0,0,alphaf;0,0,0,0;0,alphas,0,0].

If the coupling conditions are of Dirichlet and Neumann type (e.g. continuity of velocity and stress) then the preconditioned system will have two Robin conditions instead.

Parameters
matrixthe input matrix that will hold the robin coupling;
alphafparameter multiplying the b.c. of the first block;
alphasparameter multiplying the b.c. of the second block;
couplingflag specifying if we want to neglect part of the coupling (used for preconditioning purposes);
mapglobal map of the whole problem;
FESpace1FESpace of the first block;
offset1offset of the first block in the global matrix;
FESpace2FESpace of the second block;
offset2offset of the second block in the global matrix;
locDofMapstd::map with the correspondance between the numeration of the interface in the 2 FE spaces.
numerationInterfacevector containing the correspondence of the Lagrange multipliers with the interface dofs

Definition at line 167 of file MonolithicBlock.cpp.

◆ addToBlock()

void addToBlock ( const matrixPtr_Type Mat,
UInt  position 
)
virtual

adds a new block

Parameters
Matblock matrix to push
positionposition of the matrix to which we want to add Mat

Definition at line 227 of file MonolithicBlock.cpp.

◆ push_back_oper()

void push_back_oper ( MonolithicBlock Oper)
virtual

Pushes a new block.

Pushes a matrix, BCHandler, FESpace and offset in the correspondent vector

Definition at line 233 of file MonolithicBlock.cpp.

◆ blockVector()

const std::vector<matrixPtr_Type>& blockVector ( )
inline

returns the vector of pointers to the blocks (by const reference).

Definition at line 451 of file MonolithicBlock.hpp.

◆ BChVector()

const std::vector<bchandlerPtr_Type>& BChVector ( )
inline

returns the vector of pointers to the BCHandlers (by const reference).

Definition at line 457 of file MonolithicBlock.hpp.

◆ FESpaceVector()

const std::vector<fespacePtr_Type>& FESpaceVector ( )
inline

returns the vector of pointers to the FE spaces (by const reference).

Definition at line 463 of file MonolithicBlock.hpp.

◆ offsetVector()

const std::vector<UInt>& offsetVector ( )
inline

returns the vector of the offsets (by const reference).

Definition at line 469 of file MonolithicBlock.hpp.

◆ whereIsBlock()

virtual const UInt whereIsBlock ( UInt  position) const
pure virtual

◆ blockAssembling() [2/2]

virtual void blockAssembling ( const UInt  )
inlineprotectedvirtual

sums the coupling matrix in the specified position with the corresponding block

Everything (but the boundary conditions) must have been set before calling this

Reimplemented in MonolithicBlockComposed.

Definition at line 487 of file MonolithicBlock.hpp.

◆ swap()

void swap ( std::shared_ptr< Operator > &  operFrom,
std::shared_ptr< Operator > &  OperTo 
)
protected

swaps two std::shared_ptr. The tamplate argument of the shared_ptr is templated

Parameters
operFromshared_ptr to be swapped
operToshared_ptr to be swapped

Definition at line 539 of file MonolithicBlock.hpp.

◆ insert()

void insert ( std::vector< Operator > &  operFrom,
std::vector< Operator > &  OperTo 
)
protected

swaps two std::shared_ptr. The tamplate argument of the shared_ptr is templated

Parameters
operFromshared_ptr to be swapped
operToshared_ptr to be swapped

Definition at line 549 of file MonolithicBlock.hpp.

Field Documentation

◆ M_bch

std::vector<bchandlerPtr_Type> M_bch
protected

Definition at line 511 of file MonolithicBlock.hpp.

◆ M_blocks

std::vector<matrixPtr_Type> M_blocks
protected

Definition at line 512 of file MonolithicBlock.hpp.

◆ M_FESpace

std::vector<fespacePtr_Type> M_FESpace
protected

Definition at line 513 of file MonolithicBlock.hpp.

◆ M_offset

std::vector<UInt> M_offset
protected

Definition at line 514 of file MonolithicBlock.hpp.

◆ M_numerationInterface

vectorPtr_Type M_numerationInterface
protected

Definition at line 515 of file MonolithicBlock.hpp.

◆ M_comm

std::shared_ptr<Epetra_Comm> M_comm
protected

Definition at line 516 of file MonolithicBlock.hpp.


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