LifeV
FSIApplyOperator.hpp
Go to the documentation of this file.
1 /*
2  * BlockOperator.hpp
3  *
4  * Created on: Sep 20, 2010
5  * Author: uvilla
6  */
7 
8 #ifndef FSIAPPLYOPERATOR_HPP_
9 #define FSIAPPLYOPERATOR_HPP_
10 
11 #include <Epetra_Import.h>
12 #include <boost/numeric/ublas/matrix.hpp>
13 
14 #include <lifev/core/linear_algebra/BlockEpetra_Map.hpp>
15 #include <lifev/core/linear_algebra/BlockEpetra_MultiVector.hpp>
16 #include <lifev/core/linear_algebra/LinearOperatorAlgebra.hpp>
17 #include <lifev/core/array/MapEpetra.hpp>
18 #include <lifev/core/array/VectorEpetra.hpp>
19 
20 namespace LifeV
21 {
22 
23 namespace Operators
24 {
25 //! @class FSIApplyOperator
26 /*! @brief A abstract class for handling n-by-m block operators
27  * This class inherits from LifeV::LinearOperator.
28  *
29  * The Transpose is not supported yet.
30  */
31 
33 {
34 public:
35 
36  //! @name Public Typedefs
37  //@{
39  typedef super::comm_Type comm_Type;
40  typedef super::commPtr_Type commPtr_Type;
41  typedef super::map_Type map_Type;
42  typedef super::mapPtr_Type mapPtr_Type;
43  typedef super::operator_Type operator_Type;
44  typedef super::operatorPtr_Type operatorPtr_Type;
45  typedef super::vector_Type vector_Type;
46  typedef super::vectorPtr_Type vectorPtr_Type;
47 
52 
56 
57  //@}
58 
59  //! Empty Constructor
61 
62  //! @name Set Methods
63  //@{
64  //! SetUp for a "square operator"
65  /*!
66  * @param domainMap: the map of a vector in the domain of \e this
67  * is obtained by concatenating the block maps in
68  * domainMap.
69  * rangeMap is assumed to be the same of domainMap.
70  * @param comm: the communicator.
71  */
72  void setUp (const std::shared_ptr<BlockEpetra_Map> & map, const commPtr_Type & comm);
73 
74  //! SetUp for a "rectangular operator"
75  /*!
76  * @param domainMap: the map of a vector in the domain of \e this
77  * is obtained by concatenating the block maps in
78  * domainMap.
79  * @param rangeMap: the map of a vector in the range of \e this
80  * is obtained by concatenating the block maps in
81  * rangeMap.
82  * @param comm: the communicator.
83  */
84  void setUp (const std::shared_ptr<BlockEpetra_Map> & domainMap,
85  const std::shared_ptr<BlockEpetra_Map> & rangeMap,
86  const commPtr_Type & comm);
87 
88  //! SetUp when the operator is given like a boost::matrix
89  /*!
90  * @param blockOper: a dense matrix to describe the block operator
91  * @param comm: the communicator
92  */
93  void setUp (const operatorPtrContainer_Type & blockOper, const commPtr_Type & comm);
94 
95  //! set a component of the block operator
96  /*!
97  * @param iblock, jblock: The position of the block is (iblock, jblock).
98  * @param operBlock : an operator_ptr representing the block
99  */
100  void setBlock (UInt iblock, UInt jblock, const operatorPtr_Type & operBlock);
101 
102  //! Complete the block matrix with null operators
103  void fillComplete ();
104 
105  //! If true the transpose of the operator will be computed.
106  /*
107  * Not Supported yet
108  */
109  int SetUseTranspose(bool useTranspose);
110 
111  //! Set the monolithic map
112  void setMonolithicMap(const mapEpetraPtr_Type& monolithicMap){ M_monolithicMap = monolithicMap; };
113 
114  //! Compute Y = Op*X;
115  virtual int Apply(const vector_Type & X, vector_Type & Y) const;
116  //! Compute Y = Op\X;
117  /*!
118  * ApplyInverse is implemented for the matrices with block diagonal, lowerTriangular, upperTriangular form
119  */
120  virtual int ApplyInverse(const vector_Type & X, vector_Type & Y) const;
121 
122  //! Compute the Inf norm of the operator
123  /*!
124  * Not implemented yet.
125  */
126  double NormInf() const {return -1;}
127 
128  //! Returns a character string describing the operator
129  virtual const char * Label() const {return M_name.c_str();}
130 
131  //! Returns the current UseTranspose setting.
132  bool UseTranspose() const {return M_useTranspose;}
133 
134  //! Returns true if the \e this object can provide an approximate Inf-norm, false otherwise.
135  bool HasNormInf() const {return false;}
136 
137  //! Returns a pointer to the Epetra_Comm communicator associated with this operator.
138  const comm_Type & Comm() const {return *M_comm;}
139 
140  //! Returns a const pointer to the (i,j) block
141  const operatorPtr_Type& block (UInt iblock, UInt jblock) const;
142 
143  //! Returns the Epetra_Map object associated with the domain of this operator.
144  const map_Type & OperatorDomainMap() const {return *(M_domainMap->monolithicMap());}
145  //! Returns the Epetra_Map object associated with the domain of this operator as a pointer
146  const mapPtr_Type & OperatorDomainMap_ptr() const {return M_domainMap->monolithicMap();}
147  const std::shared_ptr<BlockEpetra_Map> & OperatorDomainBlockMapPtr() const {return M_domainMap;}
148 
149  //! Returns the Epetra_Map object associated with the range of this operator.
150  const map_Type & OperatorRangeMap() const {return *(M_rangeMap->monolithicMap());}
151  //! Returns the Epetra_Map object associated with the range of this operator as a pointer
152  const mapPtr_Type & OperatorRangeMap_ptr() const {return M_rangeMap->monolithicMap();}
153  const std::shared_ptr<BlockEpetra_Map> & OperatorRangeBlockMapPtr() const {return M_rangeMap;}
154 
155 protected:
156  //! Compute Y = Op*X;
157  int applyNoTranspose(const vector_Type & X, vector_Type & Y) const;
158  //! Compute Y = Op'*X;
159  int applyTranspose(const vector_Type & X, vector_Type & Y) const;
160  //! Y = diag(block(i,i)^-1)*X
161  int blockJacobi(const vector_Type & X, vector_Type & Y) const;
162  //! Y = backwardsubstitution(X)
163  int blockUpperTriangularSolve(const vector_Type & X, vector_Type & Y) const;
164  //! Y = forwardsubstitution(X)
165  int blockLowerTriangularSolve(const vector_Type & X, vector_Type & Y) const;
166  //! Change the name of the operator, (available for derivate classes).
167  void setName(const std::string & name){M_name =name;}
168 private:
169 
170  //! Number of blocks in each row
172  //! Number of blocks in each column
174 
175  //! Name of the object
176  std::string M_name;
177  //! Communicator
178  commPtr_Type M_comm;
179 
180  //! @name Maps
181  //@{
182  //! Domain Map
184  //! Range Map
186 
188  //@}
189 
190  //! block operator represented like a dense matrix of pointers to Operators
192 
193  //! whenever transpose should be used
195 };
196 
197 } /*end namespace Operators*/
198 
199 } /*end namespace */
200 #endif /* FSIAPPLYOPERATOR_HPP_ */
std::shared_ptr< mapEpetra_Type > mapEpetraPtr_Type
virtual int ApplyInverse(const vector_Type &X, vector_Type &Y) const
Compute Y = Op;.
VectorEpetra - The Epetra Vector format Wrapper.
int blockJacobi(const vector_Type &X, vector_Type &Y) const
Y = diag(block(i,i)^-1)*X.
std::vector< vectorPtr_Type > vectorPtrContainer_Type
bool M_useTranspose
whenever transpose should be used
const operatorPtr_Type & block(UInt iblock, UInt jblock) const
Returns a const pointer to the (i,j) block.
double NormInf() const
Compute the Inf norm of the operator.
const map_Type & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
commPtr_Type M_comm
Communicator.
const comm_Type & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
int blockUpperTriangularSolve(const vector_Type &X, vector_Type &Y) const
Y = backwardsubstitution(X)
Abstract class which defines the interface of a Linear Operator.
const mapPtr_Type & OperatorRangeMap_ptr() const
Returns the Epetra_Map object associated with the range of this operator as a pointer.
virtual const char * Label() const
Returns a character string describing the operator.
std::shared_ptr< BlockEpetra_Map > M_domainMap
Domain Map.
void updateInverseJacobian(const UInt &iQuadPt)
int applyNoTranspose(const vector_Type &X, vector_Type &Y) const
Compute Y = Op*X;.
bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
const std::shared_ptr< BlockEpetra_Map > & OperatorDomainBlockMapPtr() const
std::string M_name
Name of the object.
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
boost::numeric::ublas::matrix< operatorPtr_Type > operatorPtrContainer_Type
UInt M_nBlockCols
Number of blocks in each column.
int applyTranspose(const vector_Type &X, vector_Type &Y) const
Compute Y = Op&#39;*X;.
operatorPtrContainer_Type M_oper
block operator represented like a dense matrix of pointers to Operators
void setMonolithicMap(const mapEpetraPtr_Type &monolithicMap)
Set the monolithic map.
std::vector< mapPtr_Type > mapPtrContainer_Type
A abstract class for handling n-by-m block operators This class inherits from LifeV::LinearOperator.
void setUp(const operatorPtrContainer_Type &blockOper, const commPtr_Type &comm)
SetUp when the operator is given like a boost::matrix.
const std::shared_ptr< BlockEpetra_Map > & OperatorRangeBlockMapPtr() const
const map_Type & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
void setUp(const std::shared_ptr< BlockEpetra_Map > &domainMap, const std::shared_ptr< BlockEpetra_Map > &rangeMap, const commPtr_Type &comm)
SetUp for a "rectangular operator".
bool UseTranspose() const
Returns the current UseTranspose setting.
std::shared_ptr< VectorEpetra_Type > VectorEpetraPtr_Type
int blockLowerTriangularSolve(const vector_Type &X, vector_Type &Y) const
Y = forwardsubstitution(X)
void setBlock(UInt iblock, UInt jblock, const operatorPtr_Type &operBlock)
set a component of the block operator
int SetUseTranspose(bool useTranspose)
If true the transpose of the operator will be computed.
const mapPtr_Type & OperatorDomainMap_ptr() const
Returns the Epetra_Map object associated with the domain of this operator as a pointer.
UInt M_nBlockRows
Number of blocks in each row.
std::shared_ptr< vector_Type > vectorPtr_Type
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
std::shared_ptr< operator_Type > operatorPtr_Type
void fillComplete()
Complete the block matrix with null operators.
void setName(const std::string &name)
Change the name of the operator, (available for derivate classes).
virtual int Apply(const vector_Type &X, vector_Type &Y) const
Compute Y = Op*X;.
std::shared_ptr< BlockEpetra_Map > M_rangeMap
Range Map.