LifeV
LumpedOperator.hpp
Go to the documentation of this file.
1 /*
2  * LumpedOperator.hpp
3  *
4  * An Operator for Mass lumping:
5  * \hat{M}_{i,i} = M_{i,i}* \frac{\sum_{i,j} M_{i,j}}{\sum_i{M_{i,i}}
6  *
7  * Created on: Sep 17, 2010
8  * Author: uvilla
9  */
10 
11 #ifndef LUMPED_OPERATOR_HPP_
12 #define LUMPED_OPERATOR_HPP_
13 
14 #include <Epetra_CrsMatrix.h>
15 #include <Epetra_Vector.h>
16 
17 #include <lifev/core/linear_algebra/LinearOperatorAlgebra.hpp>
18 
19 namespace LifeV
20 {
21 namespace Operators
22 {
23 //! @class LumpedOperator
24 /*!
25  * @brief A class for Mass Lumping
26  * The lumped mass matrix is such that
27  * (1) The total mass is conserved
28  * (2) \hat{M}_{i,i} is proportional to M_{i,i}.
29  *
30  * The formula for the lumping is
31  * \hat{M}_{i,i} = M_{i,i}* \frac{\sum_{i,j} M_{i,j}}{\sum_i{M_{i,i}}
32  *
33  * In the case of P1 F.E. such lumping is in agreement with the
34  * one obtained using the Trapezoidal quadrature rule.
35  */
36 
38 {
39 public:
40 
42  typedef boost::shared_ptr<rowMatrix_Type> rowMatrixPtr_Type;
45 
46  //! Constructor
48  //! SetUp
49  void setUp(const rowMatrixPtr_Type & _matrix);
50  //! set the scaling constant when applying the operator
51  void setAlpha(const Real & _alpha){ M_alpha = _alpha;}
52  //! compute the lumpded operator
53  int compute();
54 
55  //! if true apply the traspose of the operator
56  int SetUseTranspose(bool UseTranspose) {M_useTranspose = UseTranspose; return 0;};
57  //! Apply the Operator
58  int Apply(const vector_Type& X, vector_Type& Y) const;
59  //! Apply the inverse of the operator
60  int ApplyInverse(const vector_Type& X, vector_Type& Y) const;
61  //! return the infinity norm of the operator
62  double NormInf() const {return M_normInf;}
63 
64  //! return the name of the operator
65  const char * Label() const{return M_name.c_str();}
66  //! return true if we are using the transpose
67  bool UseTranspose() const {return M_useTranspose;}
68  //! return true if the operator has the norm inf.
69  bool HasNormInf() const {return true;}
70  //! return the communicator
71  const comm_Type & Comm() const {return M_matrix->Comm();}
72  //! return the domain map
73  const map_Type & OperatorDomainMap() const{return M_matrix->OperatorDomainMap();}
74  //! return the range map
75  const map_Type & OperatorRangeMap() const {return M_matrix->OperatorRangeMap();}
76  //! return the lumped mass in vectorial form.
77  const lumpedMatrix_Type & getLumpedVector(){return *M_lumped;}
78  //! return the lumped mass in vectorial form (pointer).
79  const lumpedMatrixPtr_Type & getLumpedVector_ptr(){return M_lumped;}
80 
81 private:
82  //! the name of the operator
83  std::string M_name;
84  //! the matrix to be lumped
86  //! the lumped matrix in vectorial form
88  //! the infinity norm of the operator
90  //! whenever to use the transpose
92  //! the scalar coefficient for scaling the operator
94 };
95 } /*end namespace Operators*/
96 } /*end namespace */
97 #endif /* LUMPEDOPERATOR_HPP_ */
rowMatrixPtr_Type M_matrix
the matrix to be lumped
Real M_alpha
the scalar coefficient for scaling the operator
int SetUseTranspose(bool UseTranspose)
if true apply the traspose of the operator
const comm_Type & Comm() const
return the communicator
Abstract class which defines the interface of a Linear Operator.
bool M_useTranspose
whenever to use the transpose
void updateInverseJacobian(const UInt &iQuadPt)
std::string M_name
the name of the operator
boost::shared_ptr< lumpedMatrix_Type > lumpedMatrixPtr_Type
void setUp(const rowMatrixPtr_Type &_matrix)
SetUp.
const lumpedMatrix_Type & getLumpedVector()
return the lumped mass in vectorial form.
int ApplyInverse(const vector_Type &X, vector_Type &Y) const
Apply the inverse of the operator.
double NormInf() const
return the infinity norm of the operator
Real M_normInf
the infinity norm of the operator
bool HasNormInf() const
return true if the operator has the norm inf.
int Apply(const vector_Type &X, vector_Type &Y) const
Apply the Operator.
bool UseTranspose() const
return true if we are using the transpose
A class for Mass Lumping The lumped mass matrix is such that (1) The total mass is conserved (2) {M}_...
double Real
Generic real data.
Definition: LifeV.hpp:175
const map_Type & OperatorRangeMap() const
return the range map
lumpedMatrixPtr_Type M_lumped
the lumped matrix in vectorial form
boost::shared_ptr< rowMatrix_Type > rowMatrixPtr_Type
const lumpedMatrixPtr_Type & getLumpedVector_ptr()
return the lumped mass in vectorial form (pointer).
const map_Type & OperatorDomainMap() const
return the domain map
int compute()
compute the lumpded operator
const char * Label() const
return the name of the operator
void setAlpha(const Real &_alpha)
set the scaling constant when applying the operator