v0.9.0
Public Member Functions | Public Attributes | List of all members
OpAssembleMatAndVec Struct Reference
Inheritance diagram for OpAssembleMatAndVec:
[legend]
Collaboration diagram for OpAssembleMatAndVec:
[legend]

Public Member Functions

 OpAssembleMatAndVec (Mat a, Vec f)
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCoreBase::UserDataOperator
double getArea ()
 get area of face More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 get coordinates at Gauss pts. More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 get coordinates at Gauss pts (takes in account ho approx. of geometry) More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
DEPRECATED MatrixDoublegetNormalsAtGaussPt ()
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
DEPRECATED ublas::matrix_row< MatrixDoublegetNormalsAtGaussPt (const int gg)
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
const FaceElementForcesAndSourcesCoreBasegetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
DEPRECATED const FaceElementForcesAndSourcesCoreBasegetFaceElementForcesAndSourcesCore ()
 
template<int SWITCH>
MoFEMErrorCode loopSideVolumes (const string &fe_name, VolumeElementForcesAndSourcesCoreOnSideSwitch< SWITCH > &fe_method)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement ()
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

Mat A
 
Vec F
 
FTensor::Index< 'i', 3 > i
 
VectorDouble nF
 
MatrixDouble nA
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
hcurl_check_approx_in_2d.cpp.

Definition at line 92 of file hcurl_check_approx_in_2d.cpp.

Constructor & Destructor Documentation

◆ OpAssembleMatAndVec()

OpAssembleMatAndVec::OpAssembleMatAndVec ( Mat  a,
Vec  f 
)

Definition at line 95 of file hcurl_check_approx_in_2d.cpp.

96  : FaceEleOp("FIELD1", "FIELD1", OPROW | OPROWCOL), A(a), F(f) {
97  sYmm = false; // FIXME problem is symmetric, should use that
98  }
FaceEle::UserDataOperator FaceEleOp
bool sYmm
If true assume that matrix is symmetric structure.

Member Function Documentation

◆ doWork() [1/2]

MoFEMErrorCode OpAssembleMatAndVec::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
virtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 103 of file hcurl_check_approx_in_2d.cpp.

104  {
105 
107  const int nb_dofs = data.getIndices().size();
108  if (nb_dofs == 0)
110  const int nb_gauss_pts = data.getN().size1();
111  nF.resize(nb_dofs, false);
112  nF.clear();
113  auto t_base_fun = data.getFTensor1N<3>();
114  for (int gg = 0; gg != nb_gauss_pts; gg++) {
115  const double val = getArea() * getGaussPts()(2, gg);
116  for (int bb = 0; bb != nb_dofs; bb++) {
117  const double x = getCoordsAtGaussPts()(gg, 0);
118  const double y = getCoordsAtGaussPts()(gg, 1);
119  nF(bb) += val * t_base_fun(i) * ApproxFunctions::fUn(x, y)(i);
120  ++t_base_fun;
121  }
122  }
123  CHKERR VecSetValues(F, data, &*nF.data().begin(), ADD_VALUES);
125  }
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
static FTensor::Tensor1< double, 3 > fUn(const double x, const double y)
FTensor::Index< 'i', 3 > i
const VectorInt & getIndices() const
Get global indices of dofs on entity.
MoFEMErrorCode VecSetValues(Vec V, const DataForcesAndSourcesCore::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.

◆ doWork() [2/2]

MoFEMErrorCode OpAssembleMatAndVec::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)
virtual

Operator for bi-linear form, usually to calculate values on left hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 128 of file hcurl_check_approx_in_2d.cpp.

131  {
132 
134  const int nb_dofs_row = row_data.getIndices().size();
135  if (nb_dofs_row == 0)
137  const int nb_dofs_col = col_data.getIndices().size();
138  if (nb_dofs_col == 0)
140  nA.resize(nb_dofs_row, nb_dofs_col, false);
141  nA.clear();
142  const int nb_gauss_pts = row_data.getN().size1();
143  auto t_base_row = row_data.getFTensor1N<3>();
144  for (int gg = 0; gg != nb_gauss_pts; gg++) {
145  const double val = getArea() * getGaussPts()(2, gg);
146  for (int rr = 0; rr != nb_dofs_row; rr++) {
147  auto t_base_col = col_data.getFTensor1N<3>(gg, 0);
148  for (int cc = 0; cc != nb_dofs_col; cc++) {
149  nA(rr, cc) += val * t_base_row(i) * t_base_col(i);
150  ++t_base_col;
151  }
152  ++t_base_row;
153  }
154  }
155  CHKERR MatSetValues(A, row_data, col_data, &*nA.data().begin(), ADD_VALUES);
157  }
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
FTensor::Index< 'i', 3 > i
const VectorInt & getIndices() const
Get global indices of dofs on entity.
MoFEMErrorCode MatSetValues(Mat M, const DataForcesAndSourcesCore::EntData &row_data, const DataForcesAndSourcesCore::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.

Member Data Documentation

◆ A

Mat OpAssembleMatAndVec::A

Definition at line 93 of file hcurl_check_approx_in_2d.cpp.

◆ F

Vec OpAssembleMatAndVec::F

Definition at line 94 of file hcurl_check_approx_in_2d.cpp.

◆ i

FTensor::Index<'i', 3> OpAssembleMatAndVec::i

Definition at line 100 of file hcurl_check_approx_in_2d.cpp.

◆ nA

MatrixDouble OpAssembleMatAndVec::nA

Definition at line 127 of file hcurl_check_approx_in_2d.cpp.

◆ nF

VectorDouble OpAssembleMatAndVec::nF

Definition at line 102 of file hcurl_check_approx_in_2d.cpp.


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