v0.10.0
Public Member Functions | Public Attributes | List of all members
MoFEM::OpCalculateInvJacForFace Struct Reference

Calculate inverse of jacobian for face element. More...

#include <src/finite_elements/UserDataOperators.hpp>

Inheritance diagram for MoFEM::OpCalculateInvJacForFace:
[legend]
Collaboration diagram for MoFEM::OpCalculateInvJacForFace:
[legend]

Public Member Functions

 OpCalculateInvJacForFace (MatrixDouble &inv_jac)
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right 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...
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
DEPRECATED Vec getSnesF () const
 
DEPRECATED Vec getSnesX () const
 
DEPRECATED Mat getSnesA () const
 
DEPRECATED Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
double getTSaa () 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)
 
virtual ~DataOperator ()=default
 
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)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=false)
 
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

MatrixDoubleinvJac
 
- 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...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

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
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Calculate inverse of jacobian for face element.

It is assumed that face element is XY plane. Applied only for 2d problems.

Todo:
Generalize function for arbitrary face orientation in 3d space
Examples
cell_forces.cpp, contact.cpp, dynamic_elastic.cpp, eigen_elastic.cpp, hcurl_check_approx_in_2d.cpp, hcurl_divergence_operator_2d.cpp, heat_equation.cpp, helmholtz.cpp, minimal_surface_area.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, scalar_check_approximation_2d.cpp, and wave_equation.cpp.

Definition at line 1888 of file UserDataOperators.hpp.

Constructor & Destructor Documentation

◆ OpCalculateInvJacForFace()

MoFEM::OpCalculateInvJacForFace::OpCalculateInvJacForFace ( MatrixDouble inv_jac)

Definition at line 1892 of file UserDataOperators.hpp.

Member Function Documentation

◆ doWork()

MoFEMErrorCode MoFEM::OpCalculateInvJacForFace::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 149 of file UserDataOperators.cpp.

150  {
151 
153 
157 
158  size_t nb_gauss_pts = getGaussPts().size2();
159  auto &coords = getCoords();
160  double *coords_ptr = &*coords.data().begin();
161  invJac.resize(9, nb_gauss_pts, false);
162  invJac.clear();
163 
164  auto cal_inv_jac_on_tri = [&]() {
166  VectorDouble &coords = getCoords();
167  double *coords_ptr = &*coords.data().begin();
168  double j00 = 0, j01 = 0, j10 = 0, j11 = 0, j20 = 0, j21 = 0;
169 
170  // this is triangle, derivative of nodal shape functions is constant.
171  // So only need to do one node.
172  for (auto n : {0, 1, 2}) {
173  j00 += coords_ptr[3 * n + 0] * Tools::diffShapeFunMBTRI[2 * n + 0];
174  j01 += coords_ptr[3 * n + 0] * Tools::diffShapeFunMBTRI[2 * n + 1];
175  j10 += coords_ptr[3 * n + 1] * Tools::diffShapeFunMBTRI[2 * n + 0];
176  j11 += coords_ptr[3 * n + 1] * Tools::diffShapeFunMBTRI[2 * n + 1];
177 
178  // 3d
179  j20 += coords_ptr[3 * n + 2] * Tools::diffShapeFunMBTRI[2 * n + 0];
180  j21 += coords_ptr[3 * n + 2] * Tools::diffShapeFunMBTRI[2 * n + 1];
181  }
182 
183  double normal_ptr[3];
184  CHKERR Tools::getTriNormal(coords_ptr, normal_ptr);
185  const double l =
186  sqrt(normal_ptr[0] * normal_ptr[0] + normal_ptr[1] * normal_ptr[1] +
187  normal_ptr[2] * normal_ptr[2]);
188  for (auto d : {0, 1, 2})
189  normal_ptr[d] /= l;
190 
191  FTensor::Tensor2<double, 3, 3> t_jac, t_inv_jac;
192  t_jac(0, 0) = j00;
193  t_jac(0, 1) = j01;
194  t_jac(1, 0) = j10;
195  t_jac(1, 1) = j11;
196  // 3d
197  t_jac(2, 0) = j20;
198  t_jac(2, 1) = j21;
199  for (auto d : {0, 1, 2})
200  t_jac(d, d) += normal_ptr[d];
201 
202  double det;
203  CHKERR determinantTensor3by3(t_jac, det);
204  CHKERR invertTensor3by3(t_jac, det, t_inv_jac);
205 
206  auto t_inv_jac_at_pts = getFaceJac(invJac, FTensor::Number<3>());
207 
208  for (size_t gg = 0; gg != nb_gauss_pts; ++gg, ++t_inv_jac_at_pts) {
209  t_inv_jac_at_pts(i, j) = t_inv_jac(i, j);
210  }
212  };
213 
214  auto cal_inv_jac_on_quad = [&]() {
216  VectorDouble &coords = getCoords();
217  double *coords_ptr = &*coords.data().begin();
218  double *ksi_ptr = &getGaussPts()(0, 0);
219  double *zeta_ptr = &getGaussPts()(1, 0);
221 
222  auto t_inv_jac = getFaceJac(invJac, FTensor::Number<3>());
223  for (size_t gg = 0; gg != nb_gauss_pts;
224  ++gg, ++t_inv_jac, ++ksi_ptr, ++zeta_ptr) {
225  const double &ksi = *ksi_ptr;
226  const double &zeta = *zeta_ptr;
227 
228  t_jac(0, 0) = coords_ptr[3 * 0 + 0] * diffN_MBQUAD0x(zeta) +
229  coords_ptr[3 * 1 + 0] * diffN_MBQUAD1x(zeta) +
230  coords_ptr[3 * 2 + 0] * diffN_MBQUAD2x(zeta) +
231  coords_ptr[3 * 3 + 0] * diffN_MBQUAD3x(zeta);
232  t_jac(0, 1) = coords_ptr[3 * 0 + 0] * diffN_MBQUAD0y(ksi) +
233  coords_ptr[3 * 1 + 0] * diffN_MBQUAD1y(ksi) +
234  coords_ptr[3 * 2 + 0] * diffN_MBQUAD2y(ksi) +
235  coords_ptr[3 * 3 + 0] * diffN_MBQUAD3y(ksi);
236  t_jac(1, 0) = coords_ptr[3 * 0 + 1] * diffN_MBQUAD0x(zeta) +
237  coords_ptr[3 * 1 + 1] * diffN_MBQUAD1x(zeta) +
238  coords_ptr[3 * 2 + 1] * diffN_MBQUAD2x(zeta) +
239  coords_ptr[3 * 3 + 1] * diffN_MBQUAD3x(zeta);
240  t_jac(1, 1) = coords_ptr[3 * 0 + 1] * diffN_MBQUAD0y(ksi) +
241  coords_ptr[3 * 1 + 1] * diffN_MBQUAD1y(ksi) +
242  coords_ptr[3 * 2 + 1] * diffN_MBQUAD2y(ksi) +
243  coords_ptr[3 * 3 + 1] * diffN_MBQUAD3y(ksi);
244 
245  t_jac(2, 0) = coords_ptr[3 * 0 + 2] * diffN_MBQUAD0x(zeta) +
246  coords_ptr[3 * 1 + 2] * diffN_MBQUAD1x(zeta) +
247  coords_ptr[3 * 2 + 2] * diffN_MBQUAD2x(zeta) +
248  coords_ptr[3 * 3 + 2] * diffN_MBQUAD3x(zeta);
249  t_jac(2, 1) = coords_ptr[3 * 0 + 2] * diffN_MBQUAD0y(ksi) +
250  coords_ptr[3 * 1 + 2] * diffN_MBQUAD1y(ksi) +
251  coords_ptr[3 * 2 + 2] * diffN_MBQUAD2y(ksi) +
252  coords_ptr[3 * 3 + 2] * diffN_MBQUAD3y(ksi);
253 
254  FTensor::Tensor1<double, 3> t_t1{t_jac(0, 0), t_jac(1, 0), t_jac(2, 0)};
255  FTensor::Tensor1<double, 3> t_t2{t_jac(0, 1), t_jac(1, 1), t_jac(2, 1)};
256 
258  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
259  t_normal(i) /= sqrt(t_normal(j) * t_normal(j));
260  for (auto d : {0, 1, 2})
261  t_jac(d, d) += t_normal(d);
262 
263  double det;
264  CHKERR determinantTensor3by3(t_jac, det);
265  CHKERR invertTensor3by3(t_jac, det, t_inv_jac);
266 
267  }
269  };
270 
271  switch (getNumeredEntFiniteElementPtr()->getEntType()) {
272  case MBTRI:
273  CHKERR cal_inv_jac_on_tri();
274  break;
275  case MBQUAD:
276  CHKERR cal_inv_jac_on_quad();
277  break;
278  default:
279  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
280  "Operator not implemented for this entity type");
281  };
282 
284 }

Member Data Documentation

◆ invJac

MatrixDouble& MoFEM::OpCalculateInvJacForFace::invJac

Definition at line 1890 of file UserDataOperators.hpp.


The documentation for this struct was generated from the following files:
NOSPACE
@ NOSPACE
Definition: definitions.h:175
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
FTensor::Tensor2< double, 3, 3 >
MoFEM::getFaceJac
auto getFaceJac(MatrixDouble &jac, const FTensor::Number< 2 > &)
Definition: UserDataOperators.hpp:1839
diffN_MBQUAD1x
#define diffN_MBQUAD1x(y)
Definition: fem_tools.h:73
n
static Index< 'n', 3 > n
Definition: BasicFeTools.hpp:78
FTensor::d
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
diffN_MBQUAD1y
#define diffN_MBQUAD1y(x)
Definition: fem_tools.h:74
MoFEM::FaceElementForcesAndSourcesCoreBase::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:311
diffN_MBQUAD2x
#define diffN_MBQUAD2x(y)
Definition: fem_tools.h:75
MoFEM::Tools::diffShapeFunMBTRI
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:87
FTensor::Tensor1< double, 3 >
MoFEM::determinantTensor3by3
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant 3 by 3.
Definition: Templates.hpp:905
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Definition: ForcesAndSourcesCore.hpp:866
FTensor::Number
Definition: Number.hpp:12
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:604
diffN_MBQUAD3x
#define diffN_MBQUAD3x(y)
Definition: fem_tools.h:77
MoFEM::invertTensor3by3
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
Definition: Templates.hpp:886
FTensor::Index< 'i', 3 >
MoFEM::OpCalculateInvJacForFace::invJac
MatrixDouble & invJac
Definition: UserDataOperators.hpp:1890
diffN_MBQUAD2y
#define diffN_MBQUAD2y(x)
Definition: fem_tools.h:76
FTensor::levi_civita
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
Definition: Levi_Civita.hpp:617
diffN_MBQUAD0x
#define diffN_MBQUAD0x(y)
Definition: fem_tools.h:71
diffN_MBQUAD0y
#define diffN_MBQUAD0y(x)
Definition: fem_tools.h:72
diffN_MBQUAD3y
#define diffN_MBQUAD3y(x)
Definition: fem_tools.h:78
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509
i
FTensor::Index< 'i', 3 > i
Definition: matrix_function.cpp:18
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
MoFEM::Tools::getTriNormal
static MoFEMErrorCode getTriNormal(const double *coords, double *normal)
Get the Tri Normal objectGet triangle normal.
Definition: Tools.cpp:258
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:124
MoFEM::Types::VectorDouble
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:74
MoFEM::FaceElementForcesAndSourcesCoreBase::UserDataOperator::getCoords
VectorDouble & getCoords()
get triangle coordinates
Definition: FaceElementForcesAndSourcesCore.hpp:466
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Definition: ForcesAndSourcesCore.hpp:1016