v0.9.0
Public Member Functions | Public Attributes | List of all members
MagneticElement::OpStab Struct Reference

calculate and assemble stabilization matrix

\[ \mathbf{A} = \int_\Omega \epsilon \mathbf{u} \cdot \mathbf{v} \textrm{d}\Omega \]

where \(\epsilon\) is regularization parameter. More...

#include <users_modules/basic_finite_elements/magnetostatic/src/MagneticElement.hpp>

Inheritance diagram for MagneticElement::OpStab:
[legend]
Collaboration diagram for MagneticElement::OpStab:
[legend]

Public Member Functions

 OpStab (BlockData &data)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 integrate matrix More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHoGaussPtsJac ()
 
MatrixDoublegetHoGaussPtsInvJac ()
 
VectorDoublegetHoGaussPtsDetJac ()
 
auto getFTenosr0HoMeasure ()
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 Get coordinates at integration points taking geometry from field. More...
 
auto getFTensor2HoGaussPtsJac ()
 
auto getFTensor2HoGaussPtsInvJac ()
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
MoFEMErrorCode getDivergenceOfHDivBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
 Get divergence of base functions at integration point. More...
 
MoFEMErrorCode getCurlOfHCurlBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
 Get curl of base functions at integration point. More...
 
- 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 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...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
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

BlockDatablockData
 
MatrixDouble entityLocalMatrix
 
- 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

calculate and assemble stabilization matrix

\[ \mathbf{A} = \int_\Omega \epsilon \mathbf{u} \cdot \mathbf{v} \textrm{d}\Omega \]

where \(\epsilon\) is regularization parameter.

Definition at line 560 of file MagneticElement.hpp.

Constructor & Destructor Documentation

◆ OpStab()

MagneticElement::OpStab::OpStab ( BlockData data)
Examples
MagneticElement.hpp.

Definition at line 564 of file MagneticElement.hpp.

566  data.fieldName, UserDataOperator::OPROWCOL),
567  blockData(data) {
568  sYmm = true;
569  }
bool sYmm
If true assume that matrix is symmetric structure.
VolumeElementForcesAndSourcesCoreBase::UserDataOperator UserDataOperator

Member Function Documentation

◆ doWork()

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

integrate matrix

Parameters
row_sidelocal number of entity on element for row of the matrix
col_sidelocal number of entity on element for col of the matrix
row_typetype of row entity (EDGE/TRIANGLE/TETRAHEDRON)
col_typetype of col entity (EDGE/TRIANGLE/TETRAHEDRON)
row_datastructure of data, like base functions and associated methods to access those data on rows
col_datastructure of data, like base functions and associated methods to access those data on rows
Returns
error code
Examples
MagneticElement.hpp.

Definition at line 585 of file MagneticElement.hpp.

588  {
590 
591  if (row_type == MBVERTEX)
593  if (col_type == MBVERTEX)
595 
596  const int nb_row_dofs = row_data.getN().size2() / 3;
597  if (nb_row_dofs == 0)
599  const int nb_col_dofs = col_data.getN().size2() / 3;
600  if (nb_col_dofs == 0)
602  entityLocalMatrix.resize(nb_row_dofs, nb_col_dofs, false);
603  entityLocalMatrix.clear();
604 
605  if (nb_row_dofs != static_cast<int>(row_data.getFieldData().size())) {
606  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
607  "Number of base functions and DOFs on entity is different on "
608  "rows %d!=%d",
609  nb_row_dofs, row_data.getFieldData().size());
610  }
611  if (nb_col_dofs != static_cast<int>(col_data.getFieldData().size())) {
612  SETERRQ2(
613  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
614  "Number of base functions and DOFs on entity is different on cols",
615  nb_col_dofs, col_data.getFieldData().size());
616  }
617 
618  MatrixDouble row_curl_mat, col_curl_mat;
620 
621  const double c0 = 1. / blockData.mU;
622  const double c1 = blockData.ePsilon * c0;
623  const int nb_gauss_pts = row_data.getN().size1();
624 
625  for (int gg = 0; gg != nb_gauss_pts; gg++) {
626 
627  // get integration weight scaled by volume
628  double w = getGaussPts()(3, gg) * getVolume();
629  // if ho geometry is given
630  w *= getHoGaussPtsDetJac()(gg);
631 
633  &row_data.getVectorN<3>(gg)(0, HVEC0),
634  &row_data.getVectorN<3>(gg)(0, HVEC1),
635  &row_data.getVectorN<3>(gg)(0, HVEC2), 3);
636 
637  for (int aa = 0; aa != nb_row_dofs; aa++) {
638  FTensor::Tensor0<double *> t_local_mat(&entityLocalMatrix(aa, 0), 1);
640  &col_data.getVectorN<3>(gg)(0, HVEC0),
641  &col_data.getVectorN<3>(gg)(0, HVEC1),
642  &col_data.getVectorN<3>(gg)(0, HVEC2), 3);
643  for (int bb = 0; bb != nb_col_dofs; bb++) {
644  t_local_mat += c1 * w * t_row_base(i) * t_col_base(i);
645  ++t_col_base;
646  ++t_local_mat;
647  }
648  ++t_row_base;
649  }
650  }
651 
652  // cerr << entityLocalMatrix << endl;
653  // cerr << endl;
654 
655  CHKERR MatSetValues(blockData.A, nb_row_dofs, &row_data.getIndices()[0],
656  nb_col_dofs, &col_data.getIndices()[0],
657  &entityLocalMatrix(0, 0), ADD_VALUES);
658 
659  if (row_side != col_side || row_type != col_type) {
661  CHKERR MatSetValues(blockData.A, nb_col_dofs, &col_data.getIndices()[0],
662  nb_row_dofs, &row_data.getIndices()[0],
663  &entityLocalMatrix(0, 0), ADD_VALUES);
664  }
665 
667  }
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:74
double ePsilon
regularization paramater
#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
MoFEMErrorCode MatSetValues(Mat M, const DataForcesAndSourcesCore::EntData &row_data, const DataForcesAndSourcesCore::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
double mU
magnetic constant N / A2
#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

Member Data Documentation

◆ blockData

BlockData& MagneticElement::OpStab::blockData
Examples
MagneticElement.hpp.

Definition at line 563 of file MagneticElement.hpp.

◆ entityLocalMatrix

MatrixDouble MagneticElement::OpStab::entityLocalMatrix
Examples
MagneticElement.hpp.

Definition at line 571 of file MagneticElement.hpp.


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