v0.9.1
Public Member Functions | Public Attributes | List of all members
MoFEM::OpSetCovariantPiolaTransformOnFace Struct Reference

transform Hcurl base fluxes from reference element to physical triangle More...

#include <src/finite_elements/DataOperators.hpp>

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

Public Member Functions

 OpSetCovariantPiolaTransformOnFace (const VectorDouble &normal, const MatrixDouble &normals_at_pts, const VectorDouble &tangent0, const MatrixDouble &tangent0_at_pts, const VectorDouble &tangent1, const MatrixDouble &tangent1_at_pts)
 
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::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

const VectorDoublenOrmal
 
const MatrixDoublenormalsAtGaussPts
 
const VectorDoubletAngent0
 
const MatrixDoubletangent0AtGaussPt
 
const VectorDoubletAngent1
 
const MatrixDoubletangent1AtGaussPt
 
- 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...
 

Detailed Description

transform Hcurl base fluxes from reference element to physical triangle

Definition at line 522 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ OpSetCovariantPiolaTransformOnFace()

MoFEM::OpSetCovariantPiolaTransformOnFace::OpSetCovariantPiolaTransformOnFace ( const VectorDouble normal,
const MatrixDouble normals_at_pts,
const VectorDouble tangent0,
const MatrixDouble tangent0_at_pts,
const VectorDouble tangent1,
const MatrixDouble tangent1_at_pts 
)

Definition at line 531 of file DataOperators.hpp.

537  : nOrmal(normal), normalsAtGaussPts(normals_at_pts), tAngent0(tangent0),
538  tangent0AtGaussPt(tangent0_at_pts), tAngent1(tangent1),
539  tangent1AtGaussPt(tangent1_at_pts) {}
FTensor::Tensor1< double, 2 > normal(FTensor::Tensor1< T1, 3 > &t_coords, FTensor::Tensor1< T2, 2 > &t_disp)
Definition: ContactOps.hpp:128

Member Function Documentation

◆ doWork()

MoFEMErrorCode MoFEM::OpSetCovariantPiolaTransformOnFace::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 944 of file DataOperators.cpp.

945  {
947 
948  if (type != MBEDGE && type != MBTRI)
950 
954 
956  &tAngent0[0], &tAngent1[0], &nOrmal[0],
957 
958  &tAngent0[1], &tAngent1[1], &nOrmal[1],
959 
960  &tAngent0[2], &tAngent1[2], &nOrmal[2],
961 
962  3);
963  double det;
965  CHKERR determinantTensor3by3(t_m, det);
966  CHKERR invertTensor3by3(t_m, det, t_inv_m);
967 
968  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; ++b) {
969 
970  FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
971 
972  auto &baseN = data.getN(base);
973  auto &diffBaseN = data.getDiffN(base);
974 
975  int nb_dofs = baseN.size2() / 3;
976  int nb_gauss_pts = baseN.size1();
977 
978  MatrixDouble piola_n(baseN.size1(), baseN.size2());
979  MatrixDouble diff_piola_n(diffBaseN.size1(), diffBaseN.size2());
980 
981  if (nb_dofs > 0 && nb_gauss_pts > 0) {
982 
984  &baseN(0, HVEC0), &baseN(0, HVEC1), &baseN(0, HVEC2));
986  &diffBaseN(0, HVEC0_0), &diffBaseN(0, HVEC0_1),
987  &diffBaseN(0, HVEC1_0), &diffBaseN(0, HVEC1_1),
988  &diffBaseN(0, HVEC2_0), &diffBaseN(0, HVEC2_1));
989  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_transformed_h_curl(
990  &piola_n(0, HVEC0), &piola_n(0, HVEC1), &piola_n(0, HVEC2));
992  t_transformed_diff_h_curl(
993  &diff_piola_n(0, HVEC0_0), &diff_piola_n(0, HVEC0_1),
994  &diff_piola_n(0, HVEC1_0), &diff_piola_n(0, HVEC1_1),
995  &diff_piola_n(0, HVEC2_0), &diff_piola_n(0, HVEC2_1));
996 
997  int cc = 0;
998  if (normalsAtGaussPts.size1() == (unsigned int)nb_gauss_pts) {
999  // HO geometry is set, so jacobian is different at each gauss point
1001  &tangent0AtGaussPt(0, 0), &tangent1AtGaussPt(0, 0),
1002  &normalsAtGaussPts(0, 0), &tangent0AtGaussPt(0, 1),
1003  &tangent1AtGaussPt(0, 1), &normalsAtGaussPts(0, 1),
1004  &tangent0AtGaussPt(0, 2), &tangent1AtGaussPt(0, 2),
1005  &normalsAtGaussPts(0, 2), 3);
1006  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1007  CHKERR determinantTensor3by3(t_m_at_pts, det);
1008  CHKERR invertTensor3by3(t_m_at_pts, det, t_inv_m);
1009  for (int ll = 0; ll != nb_dofs; ll++) {
1010  t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
1011  t_transformed_diff_h_curl(i, k) =
1012  t_inv_m(j, i) * t_diff_h_curl(j, k);
1013  ++t_h_curl;
1014  ++t_transformed_h_curl;
1015  ++t_diff_h_curl;
1016  ++t_transformed_diff_h_curl;
1017  ++cc;
1018  }
1019  ++t_m_at_pts;
1020  }
1021  } else {
1022  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1023  for (int ll = 0; ll != nb_dofs; ll++) {
1024  t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
1025  t_transformed_diff_h_curl(i, k) =
1026  t_inv_m(j, i) * t_diff_h_curl(j, k);
1027  ++t_h_curl;
1028  ++t_transformed_h_curl;
1029  ++t_diff_h_curl;
1030  ++t_transformed_diff_h_curl;
1031  ++cc;
1032  }
1033  }
1034  }
1035  if (cc != nb_gauss_pts * nb_dofs)
1036  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "Data inconsistency");
1037 
1038  baseN.data().swap(piola_n.data());
1039  diffBaseN.data().swap(diff_piola_n.data());
1040  }
1041  }
1042 
1044 }
FTensor::Index< 'i', 3 > i
FTensor::Index< 'j', 3 > j
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant 3 by 3.
Definition: Templates.hpp:482
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:75
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:482
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513
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:463
FieldApproximationBase
approximation base
Definition: definitions.h:149
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:151
#define CHKERR
Inline error check.
Definition: definitions.h:601
FTensor::Index< 'k', 2 > k
Definition: ContactOps.hpp:28
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:412

Member Data Documentation

◆ nOrmal

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::nOrmal

Definition at line 524 of file DataOperators.hpp.

◆ normalsAtGaussPts

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::normalsAtGaussPts

Definition at line 525 of file DataOperators.hpp.

◆ tAngent0

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tAngent0

Definition at line 526 of file DataOperators.hpp.

◆ tangent0AtGaussPt

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tangent0AtGaussPt

Definition at line 527 of file DataOperators.hpp.

◆ tAngent1

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tAngent1

Definition at line 528 of file DataOperators.hpp.

◆ tangent1AtGaussPt

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tangent1AtGaussPt

Definition at line 529 of file DataOperators.hpp.


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