v0.14.0
Public Member Functions | Public Attributes | List of all members
ThermalElement::OpConvectionRhs Struct Reference

operator to calculate convection therms on body surface and assemble to rhs of equations More...

#include <users_modules/basic_finite_elements/src/ThermalElement.hpp>

Inheritance diagram for ThermalElement::OpConvectionRhs:
[legend]
Collaboration diagram for ThermalElement::OpConvectionRhs:
[legend]

Public Member Functions

 OpConvectionRhs (const std::string field_name, ConvectionData &data, CommonData &common_data, bool ho_geometry=false)
 
 OpConvectionRhs (const std::string field_name, Vec _F, ConvectionData &data, CommonData &common_data, bool ho_geometry=false)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face 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...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
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...
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, 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...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (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 () const
 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...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
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
 
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 getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. 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 doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
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

CommonDatacommonData
 
ConvectionDatadAta
 
bool hoGeometry
 
bool useTsF
 
Vec F
 
VectorDouble Nf
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- 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::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPSPACE = 1 << 3,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
- 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)>
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

operator to calculate convection therms on body surface and assemble to rhs of equations

Definition at line 492 of file ThermalElement.hpp.

Constructor & Destructor Documentation

◆ OpConvectionRhs() [1/2]

ThermalElement::OpConvectionRhs::OpConvectionRhs ( const std::string  field_name,
ConvectionData data,
CommonData common_data,
bool  ho_geometry = false 
)
inline

Definition at line 500 of file ThermalElement.hpp.

504  commonData(common_data), dAta(data), hoGeometry(ho_geometry),
505  useTsF(true) {}

◆ OpConvectionRhs() [2/2]

ThermalElement::OpConvectionRhs::OpConvectionRhs ( const std::string  field_name,
Vec  _F,
ConvectionData data,
CommonData common_data,
bool  ho_geometry = false 
)
inline

Definition at line 508 of file ThermalElement.hpp.

512  commonData(common_data), dAta(data), hoGeometry(ho_geometry),
513  useTsF(false), F(_F) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode ThermalElement::OpConvectionRhs::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)

brief calculate Convection condition on the right hand side R=int_S N^T*alpha*N_f dS

Definition at line 345 of file ThermalElement.cpp.

346  {
348 
349  if (dAta.tRis.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
350  dAta.tRis.end()) {
352  }
353  if (data.getIndices().size() == 0)
355 
356  const auto &dof_ptr = data.getFieldDofs()[0];
357  int rank = dof_ptr->getNbOfCoeffs();
358 
359  int nb_row_dofs = data.getIndices().size() / rank;
360 
361  Nf.resize(data.getIndices().size());
362  Nf.clear();
363 
364  for (unsigned int gg = 0; gg < data.getN().size1(); gg++) {
365 
366  double T_at_Gauss_pt = commonData.temperatureAtGaussPts[gg];
367  double convectionConst;
368  if (hoGeometry) {
369  double area = norm_2(getNormalsAtGaussPts(gg)) * 0.5;
370  convectionConst =
371  dAta.cOnvection * area * (T_at_Gauss_pt - dAta.tEmperature);
372  } else {
373  convectionConst =
374  dAta.cOnvection * getArea() * (T_at_Gauss_pt - dAta.tEmperature);
375  }
376  double val = getGaussPts()(2, gg) * convectionConst;
377  ublas::noalias(Nf) += val * data.getN(gg, nb_row_dofs);
378  }
379 
380  if (useTsF) {
381  CHKERR VecSetValues(getFEMethod()->ts_F, data.getIndices().size(),
382  &data.getIndices()[0], &Nf[0], ADD_VALUES);
383  } else {
384  CHKERR VecSetValues(F, data.getIndices().size(), &data.getIndices()[0],
385  &Nf[0], ADD_VALUES);
386  }
387 
389 }

Member Data Documentation

◆ commonData

CommonData& ThermalElement::OpConvectionRhs::commonData

Definition at line 496 of file ThermalElement.hpp.

◆ dAta

ConvectionData& ThermalElement::OpConvectionRhs::dAta

Definition at line 497 of file ThermalElement.hpp.

◆ F

Vec ThermalElement::OpConvectionRhs::F

Definition at line 507 of file ThermalElement.hpp.

◆ hoGeometry

bool ThermalElement::OpConvectionRhs::hoGeometry

Definition at line 498 of file ThermalElement.hpp.

◆ Nf

VectorDouble ThermalElement::OpConvectionRhs::Nf

Definition at line 515 of file ThermalElement.hpp.

◆ useTsF

bool ThermalElement::OpConvectionRhs::useTsF

Definition at line 499 of file ThermalElement.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
ThermalElement::ConvectionData::tRis
Range tRis
Definition: ThermalElement.hpp:120
ThermalElement::OpConvectionRhs::F
Vec F
Definition: ThermalElement.hpp:507
ThermalElement::OpConvectionRhs::Nf
VectorDouble Nf
Definition: ThermalElement.hpp:515
_F
#define _F(n)
Definition: quad.c:25
ThermalElement::OpConvectionRhs::dAta
ConvectionData & dAta
Definition: ThermalElement.hpp:497
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1589
MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Definition: ForcesAndSourcesCore.hpp:1236
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPts
MatrixDouble & getNormalsAtGaussPts()
if higher order geometry return normals at Gauss pts.
Definition: FaceElementForcesAndSourcesCore.hpp:290
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:86
ThermalElement::ConvectionData::tEmperature
double tEmperature
Definition: ThermalElement.hpp:118
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
ThermalElement::ConvectionData::cOnvection
double cOnvection
Definition: ThermalElement.hpp:117
ThermalElement::OpConvectionRhs::commonData
CommonData & commonData
Definition: ThermalElement.hpp:496
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Definition: ForcesAndSourcesCore.hpp:1000
MoFEM::EntitiesFieldData::EntData::getFieldDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
Definition: EntitiesFieldData.hpp:1269
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Definition: ForcesAndSourcesCore.hpp:1042
ThermalElement::CommonData::temperatureAtGaussPts
VectorDouble temperatureAtGaussPts
Definition: ThermalElement.hpp:145
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
ThermalElement::OpConvectionRhs::hoGeometry
bool hoGeometry
Definition: ThermalElement.hpp:498
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getArea
double getArea()
get area of face
Definition: FaceElementForcesAndSourcesCore.hpp:239
ThermalElement::OpConvectionRhs::useTsF
bool useTsF
Definition: ThermalElement.hpp:499
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
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567