v0.14.0
Public Member Functions | Public Attributes | List of all members
ThermalElement::OpGetFieldAtGaussPts< OP > Struct Template Reference

operator to calculate temperature and rate of temperature at Gauss points More...

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

Inheritance diagram for ThermalElement::OpGetFieldAtGaussPts< OP >:
[legend]
Collaboration diagram for ThermalElement::OpGetFieldAtGaussPts< OP >:
[legend]

Public Member Functions

 OpGetFieldAtGaussPts (const std::string field_name, VectorDouble &field_at_gauss_pts)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 operator calculating temperature and rate of temperature More...
 

Public Attributes

VectorDouble & fieldAtGaussPts
 

Detailed Description

template<typename OP>
struct ThermalElement::OpGetFieldAtGaussPts< OP >

operator to calculate temperature and rate of temperature at Gauss points

Definition at line 177 of file ThermalElement.hpp.

Constructor & Destructor Documentation

◆ OpGetFieldAtGaussPts()

template<typename OP >
ThermalElement::OpGetFieldAtGaussPts< OP >::OpGetFieldAtGaussPts ( const std::string  field_name,
VectorDouble &  field_at_gauss_pts 
)
inline

Definition at line 180 of file ThermalElement.hpp.

182  : OP::UserDataOperator(field_name, OP::UserDataOperator::OPROW),
183  fieldAtGaussPts(field_at_gauss_pts) {}

Member Function Documentation

◆ doWork()

template<typename OP >
MoFEMErrorCode ThermalElement::OpGetFieldAtGaussPts< OP >::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
inline

operator calculating temperature and rate of temperature

temperature temperature or rate of temperature is calculated multiplying shape functions by degrees of freedom

Definition at line 190 of file ThermalElement.hpp.

191  {
193  try {
194 
195  if (data.getFieldData().size() == 0)
197  int nb_dofs = data.getFieldData().size();
198  int nb_gauss_pts = data.getN().size1();
199 
200  // initialize
201  fieldAtGaussPts.resize(nb_gauss_pts);
202  if (type == MBVERTEX) {
203  // loop over shape functions on entities always start from
204  // vertices, so if nodal shape functions are processed, vector of
205  // field values is zero at initialization
206  std::fill(fieldAtGaussPts.begin(), fieldAtGaussPts.end(), 0);
207  }
208 
209  for (int gg = 0; gg < nb_gauss_pts; gg++) {
210  fieldAtGaussPts[gg] +=
211  inner_prod(data.getN(gg, nb_dofs), data.getFieldData());
212  }
213 
214  } catch (const std::exception &ex) {
215  std::ostringstream ss;
216  ss << "throw in method: " << ex.what() << std::endl;
217  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
218  }
219 
221  }

Member Data Documentation

◆ fieldAtGaussPts

template<typename OP >
VectorDouble& ThermalElement::OpGetFieldAtGaussPts< OP >::fieldAtGaussPts

Definition at line 179 of file ThermalElement.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MOFEM_STD_EXCEPTION_THROW
@ MOFEM_STD_EXCEPTION_THROW
Definition: definitions.h:39
convert.type
type
Definition: convert.py:64
ThermalElement::OpGetFieldAtGaussPts::fieldAtGaussPts
VectorDouble & fieldAtGaussPts
Definition: ThermalElement.hpp:179
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440