v0.14.0
Public Member Functions | Private Attributes | List of all members
OpEnergyDensity Struct Reference

#include <tutorials/scl-12/src/electrostatics.hpp>

Inheritance diagram for OpEnergyDensity:
[legend]
Collaboration diagram for OpEnergyDensity:
[legend]

Public Member Functions

 OpEnergyDensity (const std::string &field_name, boost::shared_ptr< MatrixDouble > grad_u_negative, boost::shared_ptr< VectorDouble > energy_densiity_ptr, boost::shared_ptr< std::map< int, BlockData >> perm_block_sets_ptr, boost::shared_ptr< DataAtIntegrationPts > common_data_ptr)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 

Private Attributes

boost::shared_ptr< MatrixDouble > gradUNegative
 
boost::shared_ptr< VectorDouble > energyDensity
 
boost::shared_ptr< std::map< int, BlockData > > permBlockSetsPtr
 
boost::shared_ptr< DataAtIntegrationPtscommonDataPtr
 

Detailed Description

Definition at line 195 of file electrostatics.hpp.

Constructor & Destructor Documentation

◆ OpEnergyDensity()

OpEnergyDensity::OpEnergyDensity ( const std::string &  field_name,
boost::shared_ptr< MatrixDouble >  grad_u_negative,
boost::shared_ptr< VectorDouble >  energy_densiity_ptr,
boost::shared_ptr< std::map< int, BlockData >>  perm_block_sets_ptr,
boost::shared_ptr< DataAtIntegrationPts common_data_ptr 
)
inline

Definition at line 197 of file electrostatics.hpp.

203  : DomainEleOp(field_name, field_name, OPROWCOL, false),
204  energyDensity(energy_densiity_ptr), gradUNegative(grad_u_negative),
205  permBlockSetsPtr(perm_block_sets_ptr), commonDataPtr(common_data_ptr) {
206  std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
207  doEntities[MBVERTEX] = true;
208  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpEnergyDensity::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
inline

Definition at line 210 of file electrostatics.hpp.

213  {
215 
216  size_t nb_gauss_pts = getGaussPts().size2();
217  energyDensity->resize(nb_gauss_pts, false);
218  energyDensity->clear();
219  auto t_energy_density = getFTensor0FromVec(*energyDensity);
220  auto t_negative_grad_u = getFTensor1FromMat<SPACE_DIM>(*gradUNegative);
222  double blockPermittivity = 0.0;
223  for (const auto &n : *permBlockSetsPtr) {
224  if (n.second.domainEnts.find(getFEEntityHandle()) !=
225  n.second.domainEnts.end()) {
226  blockPermittivity = n.second.epsPermit;
227  }
228  }
229  // E = 0.5 * E^2 * Perm
230  for (int gg = 0; gg != nb_gauss_pts; gg++) {
231  t_energy_density =
232  0.5 * t_negative_grad_u(I) * t_negative_grad_u(I) * blockPermittivity;
233  ++t_negative_grad_u;
234  ++t_energy_density;
235  }
236 
238  }

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<DataAtIntegrationPts> OpEnergyDensity::commonDataPtr
private

Definition at line 244 of file electrostatics.hpp.

◆ energyDensity

boost::shared_ptr<VectorDouble> OpEnergyDensity::energyDensity
private

Definition at line 242 of file electrostatics.hpp.

◆ gradUNegative

boost::shared_ptr<MatrixDouble> OpEnergyDensity::gradUNegative
private

Definition at line 241 of file electrostatics.hpp.

◆ permBlockSetsPtr

boost::shared_ptr<std::map<int, BlockData> > OpEnergyDensity::permBlockSetsPtr
private

Definition at line 243 of file electrostatics.hpp.


The documentation for this struct was generated from the following file:
DomainEleOp
DomainEle::UserDataOperator DomainEleOp
Definition: electrostatics.hpp:16
OpEnergyDensity::gradUNegative
boost::shared_ptr< MatrixDouble > gradUNegative
Definition: electrostatics.hpp:241
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
OpEnergyDensity::energyDensity
boost::shared_ptr< VectorDouble > energyDensity
Definition: electrostatics.hpp:242
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index
Definition: Index.hpp:23
convert.n
n
Definition: convert.py:82
OpEnergyDensity::commonDataPtr
boost::shared_ptr< DataAtIntegrationPts > commonDataPtr
Definition: electrostatics.hpp:244
OpEnergyDensity::permBlockSetsPtr
boost::shared_ptr< std::map< int, BlockData > > permBlockSetsPtr
Definition: electrostatics.hpp:243
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