v0.8.20
Public Member Functions | Public Attributes | List of all members
MoFEM::OpSetCovariantPiolaTransformOnTriangle Struct Reference

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

#include <src/finite_elements/DataOperators.hpp>

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

Public Member Functions

 OpSetCovariantPiolaTransformOnTriangle (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, 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 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

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...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 

Detailed Description

transform Hcurl base fluxes from reference element to physical triangle

Definition at line 679 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ OpSetCovariantPiolaTransformOnTriangle()

MoFEM::OpSetCovariantPiolaTransformOnTriangle::OpSetCovariantPiolaTransformOnTriangle ( 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 688 of file DataOperators.hpp.

694  : nOrmal(normal), normalsAtGaussPts(normals_at_pts), tAngent0(tangent0),
695  tangent0AtGaussPt(tangent0_at_pts), tAngent1(tangent1),
696  tangent1AtGaussPt(tangent1_at_pts) {}

Member Function Documentation

◆ doWork()

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

1409  {
1410 
1412 
1413  if (type != MBEDGE && type != MBTRI)
1415 
1419 
1420  // double zero = 0;
1422  &tAngent0[0], &tAngent1[0], &nOrmal[0], &tAngent0[1], &tAngent1[1],
1423  &nOrmal[1], &tAngent0[2], &tAngent1[2], &nOrmal[2], 3);
1424  double det;
1426  ierr = determinantTensor3by3(t_m, det);
1427  CHKERRG(ierr);
1428  ierr = invertTensor3by3(t_m, det, t_inv_m);
1429  CHKERRG(ierr);
1430 
1431  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
1432 
1433  FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
1434 
1435  int nb_dofs = data.getN(base).size2() / 3;
1436  int nb_gauss_pts = data.getN(base).size1();
1437 
1438  MatrixDouble piola_n(data.getN(base).size1(), data.getN(base).size2());
1439  MatrixDouble diff_piola_n(data.getDiffN(base).size1(),
1440  data.getDiffN(base).size2());
1441 
1442  if (nb_dofs > 0 && nb_gauss_pts > 0) {
1443 
1444  FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
1446  &data.getN(base)(0, HVEC0), &data.getN(base)(0, HVEC1),
1447  &data.getN(base)(0, HVEC2));
1449  &data.getDiffN(base)(0, HVEC0_0), &data.getDiffN(base)(0, HVEC0_1),
1450  &data.getDiffN(base)(0, HVEC1_0), &data.getDiffN(base)(0, HVEC1_1),
1451  &data.getDiffN(base)(0, HVEC2_0), &data.getDiffN(base)(0, HVEC2_1));
1452  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_transformed_h_curl(
1453  &piola_n(0, HVEC0), &piola_n(0, HVEC1), &piola_n(0, HVEC2));
1455  t_transformed_diff_h_curl(
1456  &diff_piola_n(0, HVEC0_0), &diff_piola_n(0, HVEC0_1),
1457  &diff_piola_n(0, HVEC1_0), &diff_piola_n(0, HVEC1_1),
1458  &diff_piola_n(0, HVEC2_0), &diff_piola_n(0, HVEC2_1));
1459 
1460  int cc = 0;
1461  if (normalsAtGaussPts.size1() == (unsigned int)nb_gauss_pts) {
1462  // HO geometry is set, so jacobian is different at each gauss point
1464  &tangent0AtGaussPt(0, 0), &tangent1AtGaussPt(0, 0),
1465  &normalsAtGaussPts(0, 0), &tangent0AtGaussPt(0, 1),
1466  &tangent1AtGaussPt(0, 1), &normalsAtGaussPts(0, 1),
1467  &tangent0AtGaussPt(0, 2), &tangent1AtGaussPt(0, 2),
1468  &normalsAtGaussPts(0, 2), 3);
1469  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1470  ierr = determinantTensor3by3(t_m_at_pts, det);
1471  CHKERRG(ierr);
1472  ierr = invertTensor3by3(t_m_at_pts, det, t_inv_m);
1473  CHKERRG(ierr);
1474  for (int ll = 0; ll != nb_dofs; ll++) {
1475  t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
1476  t_transformed_diff_h_curl(i, k) =
1477  t_inv_m(j, i) * t_diff_h_curl(j, k);
1478  ++t_h_curl;
1479  ++t_transformed_h_curl;
1480  ++t_diff_h_curl;
1481  ++t_transformed_diff_h_curl;
1482  ++cc;
1483  }
1484  ++t_m_at_pts;
1485  }
1486  } else {
1487  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1488  for (int ll = 0; ll != nb_dofs; ll++) {
1489  t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
1490  t_transformed_diff_h_curl(i, k) =
1491  t_inv_m(j, i) * t_diff_h_curl(j, k);
1492  ++t_h_curl;
1493  ++t_transformed_h_curl;
1494  ++t_diff_h_curl;
1495  ++t_transformed_diff_h_curl;
1496  ++cc;
1497  }
1498  }
1499  }
1500  if (cc != nb_gauss_pts * nb_dofs) {
1501  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "Data inconsistency");
1502  }
1503  data.getN(base).data().swap(piola_n.data());
1504  data.getDiffN(base).data().swap(diff_piola_n.data());
1505  }
1506  }
1507 
1509 }
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.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:77
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant.
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:542
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
FieldApproximationBase
approximation base
Definition: definitions.h:141
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:143

Member Data Documentation

◆ nOrmal

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnTriangle::nOrmal

Definition at line 681 of file DataOperators.hpp.

◆ normalsAtGaussPts

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnTriangle::normalsAtGaussPts

Definition at line 682 of file DataOperators.hpp.

◆ tAngent0

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnTriangle::tAngent0

Definition at line 683 of file DataOperators.hpp.

◆ tangent0AtGaussPt

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnTriangle::tangent0AtGaussPt

Definition at line 684 of file DataOperators.hpp.

◆ tAngent1

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnTriangle::tAngent1

Definition at line 685 of file DataOperators.hpp.

◆ tangent1AtGaussPt

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnTriangle::tangent1AtGaussPt

Definition at line 686 of file DataOperators.hpp.


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