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

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

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

Public Member Functions

 OpTotalEnergy (const std::string &field_name, boost::shared_ptr< MatrixDouble > grad_u_negative, boost::shared_ptr< std::map< int, BlockData >> perm_block_sets_ptr, boost::shared_ptr< std::map< int, BlockData >> int_domain_block_sets_ptr, boost::shared_ptr< DataAtIntegrationPts > common_data_ptr, SmartPetscObj< Vec > petsc_vec_energy)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 

Private Attributes

boost::shared_ptr< MatrixDouble > gradUNegative
 
boost::shared_ptr< std::map< int, BlockData > > permBlockSetsPtr
 
boost::shared_ptr< std::map< int, BlockData > > intDomainBlocksetsPtr
 
boost::shared_ptr< DataAtIntegrationPtscommonDataPtr
 
SmartPetscObj< Vec > petscTotalEnergy
 

Detailed Description

Definition at line 135 of file electrostatics.hpp.

Constructor & Destructor Documentation

◆ OpTotalEnergy()

OpTotalEnergy::OpTotalEnergy ( const std::string &  field_name,
boost::shared_ptr< MatrixDouble >  grad_u_negative,
boost::shared_ptr< std::map< int, BlockData >>  perm_block_sets_ptr,
boost::shared_ptr< std::map< int, BlockData >>  int_domain_block_sets_ptr,
boost::shared_ptr< DataAtIntegrationPts common_data_ptr,
SmartPetscObj< Vec >  petsc_vec_energy 
)
inline

Definition at line 136 of file electrostatics.hpp.

143  : DomainEleOp(field_name, DomainEleOp::OPROW),
144  gradUNegative(grad_u_negative), permBlockSetsPtr(perm_block_sets_ptr),
145  intDomainBlocksetsPtr(int_domain_block_sets_ptr),
146  commonDataPtr(common_data_ptr), petscTotalEnergy(petsc_vec_energy) {
147  std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
148  doEntities[MBVERTEX] = true;
149  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpTotalEnergy::doWork ( int  side,
EntityType  type,
EntData data 
)
inline

Definition at line 150 of file electrostatics.hpp.

150  {
152 
153  auto t_negative_grad_u = getFTensor1FromMat<SPACE_DIM>(*gradUNegative);
155  double area = getMeasure();
156  auto t_w = getFTensor0IntegrationWeight();
157  auto nb_gauss_pts = getGaussPts().size2();
158  double totalEnergy = 0.0;
159  int index = 0;
160  // Total Energy Density = 0.5 * E^2 * Perm * Area
161  for (const auto &m : *intDomainBlocksetsPtr) {
162 
163  double blockPermittivity = 0.0;
164 
165  if (m.second.internalDomainEnts.find(getFEEntityHandle()) !=
166  m.second.internalDomainEnts.end()) {
167  for (const auto &n : *permBlockSetsPtr) {
168  if (n.second.domainEnts.find(getFEEntityHandle()) !=
169  n.second.domainEnts.end()) {
170  blockPermittivity = n.second.epsPermit;
171  }
172  }
173 
174  for (int gg = 0; gg != nb_gauss_pts; gg++) {
175  totalEnergy += 0.5 * t_negative_grad_u(I) * t_negative_grad_u(I) *
176  blockPermittivity * t_w * area;
177  ++t_negative_grad_u;
178  ++t_w;
179  }
180  }
181  }
182  CHKERR VecSetValues(petscTotalEnergy, 1, &index, &totalEnergy, ADD_VALUES);
183 
185  }

Member Data Documentation

◆ commonDataPtr

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

Definition at line 191 of file electrostatics.hpp.

◆ gradUNegative

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

Definition at line 188 of file electrostatics.hpp.

◆ intDomainBlocksetsPtr

boost::shared_ptr<std::map<int, BlockData> > OpTotalEnergy::intDomainBlocksetsPtr
private

Definition at line 190 of file electrostatics.hpp.

◆ permBlockSetsPtr

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

Definition at line 189 of file electrostatics.hpp.

◆ petscTotalEnergy

SmartPetscObj<Vec> OpTotalEnergy::petscTotalEnergy
private

Definition at line 192 of file electrostatics.hpp.


The documentation for this struct was generated from the following file:
DomainEleOp
DomainEle::UserDataOperator DomainEleOp
Definition: electrostatics.hpp:16
OpTotalEnergy::petscTotalEnergy
SmartPetscObj< Vec > petscTotalEnergy
Definition: electrostatics.hpp:192
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1589
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
OpTotalEnergy::permBlockSetsPtr
boost::shared_ptr< std::map< int, BlockData > > permBlockSetsPtr
Definition: electrostatics.hpp:189
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
OpTotalEnergy::commonDataPtr
boost::shared_ptr< DataAtIntegrationPts > commonDataPtr
Definition: electrostatics.hpp:191
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
OpTotalEnergy::gradUNegative
boost::shared_ptr< MatrixDouble > gradUNegative
Definition: electrostatics.hpp:188
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
OpTotalEnergy::intDomainBlocksetsPtr
boost::shared_ptr< std::map< int, BlockData > > intDomainBlocksetsPtr
Definition: electrostatics.hpp:190
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