v0.14.0
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, EntitiesFieldData::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, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &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
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
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::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 

Detailed Description

transform Hcurl base fluxes from reference element to physical triangle

Deprecated:
It is used in contact elements. Contact elements should be minified to work as face element,

Definition at line 457 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 
)
inline

Definition at line 466 of file DataOperators.hpp.

472  : nOrmal(normal), normalsAtGaussPts(normals_at_pts), tAngent0(tangent0),
473  tangent0AtGaussPt(tangent0_at_pts), tAngent1(tangent1),
474  tangent1AtGaussPt(tangent1_at_pts) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode MoFEM::OpSetCovariantPiolaTransformOnFace::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
virtual

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

Reimplemented from MoFEM::DataOperator.

Definition at line 590 of file DataOperators.cpp.

591  {
593 
594  const auto type_dim = moab::CN::Dimension(type);
595  if (type_dim != 1 && type_dim != 2)
597 
601 
603  &tAngent0[0], &tAngent1[0], &nOrmal[0],
604 
605  &tAngent0[1], &tAngent1[1], &nOrmal[1],
606 
607  &tAngent0[2], &tAngent1[2], &nOrmal[2]);
608  double det;
610  CHKERR determinantTensor3by3(t_m, det);
611  CHKERR invertTensor3by3(t_m, det, t_inv_m);
612 
613  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; ++b) {
614 
615  FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
616 
617  auto &baseN = data.getN(base);
618  auto &diffBaseN = data.getDiffN(base);
619 
620  int nb_dofs = baseN.size2() / 3;
621  int nb_gauss_pts = baseN.size1();
622 
623  MatrixDouble piola_n(baseN.size1(), baseN.size2());
624  MatrixDouble diff_piola_n(diffBaseN.size1(), diffBaseN.size2());
625 
626  if (nb_dofs > 0 && nb_gauss_pts > 0) {
627 
629  &baseN(0, HVEC0), &baseN(0, HVEC1), &baseN(0, HVEC2));
631  &diffBaseN(0, HVEC0_0), &diffBaseN(0, HVEC0_1),
632  &diffBaseN(0, HVEC1_0), &diffBaseN(0, HVEC1_1),
633  &diffBaseN(0, HVEC2_0), &diffBaseN(0, HVEC2_1));
634  FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_transformed_h_curl(
635  &piola_n(0, HVEC0), &piola_n(0, HVEC1), &piola_n(0, HVEC2));
637  t_transformed_diff_h_curl(
638  &diff_piola_n(0, HVEC0_0), &diff_piola_n(0, HVEC0_1),
639  &diff_piola_n(0, HVEC1_0), &diff_piola_n(0, HVEC1_1),
640  &diff_piola_n(0, HVEC2_0), &diff_piola_n(0, HVEC2_1));
641 
642  int cc = 0;
643  if (normalsAtGaussPts.size1() == (unsigned int)nb_gauss_pts) {
644  // HO geometry is set, so jacobian is different at each gauss point
646  &tangent0AtGaussPt(0, 0), &tangent1AtGaussPt(0, 0),
647  &normalsAtGaussPts(0, 0), &tangent0AtGaussPt(0, 1),
648  &tangent1AtGaussPt(0, 1), &normalsAtGaussPts(0, 1),
649  &tangent0AtGaussPt(0, 2), &tangent1AtGaussPt(0, 2),
650  &normalsAtGaussPts(0, 2));
651  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
652  CHKERR determinantTensor3by3(t_m_at_pts, det);
653  CHKERR invertTensor3by3(t_m_at_pts, det, t_inv_m);
654  for (int ll = 0; ll != nb_dofs; ll++) {
655  t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
656  t_transformed_diff_h_curl(i, k) =
657  t_inv_m(j, i) * t_diff_h_curl(j, k);
658  ++t_h_curl;
659  ++t_transformed_h_curl;
660  ++t_diff_h_curl;
661  ++t_transformed_diff_h_curl;
662  ++cc;
663  }
664  ++t_m_at_pts;
665  }
666  } else {
667  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
668  for (int ll = 0; ll != nb_dofs; ll++) {
669  t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
670  t_transformed_diff_h_curl(i, k) =
671  t_inv_m(j, i) * t_diff_h_curl(j, k);
672  ++t_h_curl;
673  ++t_transformed_h_curl;
674  ++t_diff_h_curl;
675  ++t_transformed_diff_h_curl;
676  ++cc;
677  }
678  }
679  }
680  if (cc != nb_gauss_pts * nb_dofs)
681  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "Data inconsistency");
682 
683  baseN.swap(piola_n);
684  diffBaseN.swap(diff_piola_n);
685  }
686  }
687 
689 }

Member Data Documentation

◆ nOrmal

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::nOrmal

Definition at line 459 of file DataOperators.hpp.

◆ normalsAtGaussPts

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::normalsAtGaussPts

Definition at line 460 of file DataOperators.hpp.

◆ tAngent0

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tAngent0

Definition at line 461 of file DataOperators.hpp.

◆ tangent0AtGaussPt

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tangent0AtGaussPt

Definition at line 462 of file DataOperators.hpp.

◆ tAngent1

const VectorDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tAngent1

Definition at line 463 of file DataOperators.hpp.

◆ tangent1AtGaussPt

const MatrixDouble& MoFEM::OpSetCovariantPiolaTransformOnFace::tangent1AtGaussPt

Definition at line 464 of file DataOperators.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
LASTBASE
@ LASTBASE
Definition: definitions.h:69
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::OpSetCovariantPiolaTransformOnFace::normalsAtGaussPts
const MatrixDouble & normalsAtGaussPts
Definition: DataOperators.hpp:460
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
HVEC0_1
@ HVEC0_1
Definition: definitions.h:208
MoFEM::OpSetCovariantPiolaTransformOnFace::nOrmal
const VectorDouble & nOrmal
Definition: DataOperators.hpp:459
HVEC0_0
@ HVEC0_0
Definition: definitions.h:205
HVEC1_1
@ HVEC1_1
Definition: definitions.h:209
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
HVEC1
@ HVEC1
Definition: definitions.h:199
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
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:1588
HVEC2_1
@ HVEC2_1
Definition: definitions.h:210
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
convert.type
type
Definition: convert.py:64
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::OpSetCovariantPiolaTransformOnFace::tangent1AtGaussPt
const MatrixDouble & tangent1AtGaussPt
Definition: DataOperators.hpp:464
FTensor::Index< 'i', 3 >
MoFEM::determinantTensor3by3
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1540
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
HVEC2_0
@ HVEC2_0
Definition: definitions.h:207
MoFEM::OpSetCovariantPiolaTransformOnFace::tangent0AtGaussPt
const MatrixDouble & tangent0AtGaussPt
Definition: DataOperators.hpp:462
MoFEM::OpSetCovariantPiolaTransformOnFace::tAngent0
const VectorDouble & tAngent0
Definition: DataOperators.hpp:461
MoFEM::OpSetCovariantPiolaTransformOnFace::tAngent1
const VectorDouble & tAngent1
Definition: DataOperators.hpp:463
HVEC2
@ HVEC2
Definition: definitions.h:199
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:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
HVEC1_0
@ HVEC1_0
Definition: definitions.h:206
HVEC0
@ HVEC0
Definition: definitions.h:199