v0.15.0
Loading...
Searching...
No Matches
OpAssembleG Struct Reference

Assemble G *. More...

#include "tutorials/cor-7/src/ElasticityMixedFormulation.hpp"

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

Public Member Functions

 OpAssembleG (DataAtIntegrationPts &common_data, BlockData &data)
 
PetscErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 

Public Attributes

DataAtIntegrationPtscommonData
 
MatrixDouble locG
 
BlockDatadAta
 

Detailed Description

Assemble G *.

  • \[ {\bf{G}} = - \int\limits_\Omega {{{\bf{B}}^T}{\bf m}{{\bf{N}}_p}d\Omega } \]

Examples
ElasticityMixedFormulation.hpp, and elasticity_mixed_formulation.cpp.

Definition at line 214 of file ElasticityMixedFormulation.hpp.

Constructor & Destructor Documentation

◆ OpAssembleG()

OpAssembleG::OpAssembleG ( DataAtIntegrationPts & common_data,
BlockData & data )
inline
Examples
ElasticityMixedFormulation.hpp.

Definition at line 221 of file ElasticityMixedFormulation.hpp.

222 : VolumeElementForcesAndSourcesCore::UserDataOperator(
223 "U", "P", UserDataOperator::OPROWCOL),
224 commonData(common_data), dAta(data) {
225 sYmm = false;
226 }
DataAtIntegrationPts & commonData

Member Function Documentation

◆ doWork()

PetscErrorCode OpAssembleG::doWork ( int row_side,
int col_side,
EntityType row_type,
EntityType col_type,
EntitiesFieldData::EntData & row_data,
EntitiesFieldData::EntData & col_data )
inline
Examples
ElasticityMixedFormulation.hpp.

Definition at line 228 of file ElasticityMixedFormulation.hpp.

231 {
232
234
235 const int row_nb_dofs = row_data.getIndices().size();
236 if (!row_nb_dofs)
238 const int col_nb_dofs = col_data.getIndices().size();
239 if (!col_nb_dofs)
241
242 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
243 dAta.tEts.end()) {
245 }
247
248 // Set size can clear local tangent matrix
249 locG.resize(row_nb_dofs, col_nb_dofs, false);
250 locG.clear();
251 const int row_nb_gauss_pts = row_data.getN().size1();
252 if (!row_nb_gauss_pts)
254 const int row_nb_base_functions = row_data.getN().size2();
255 auto row_diff_base_functions = row_data.getFTensor1DiffN<3>();
256
258 FTensor::Index<'i', 3> i;
259 FTensor::Index<'j', 3> j;
260
261 // INTEGRATION
262 for (int gg = 0; gg != row_nb_gauss_pts; gg++) {
263
264 // Get volume and integration weight
265 double w = getVolume() * getGaussPts()(3, gg);
266
267 int row_bb = 0;
268 for (; row_bb != row_nb_dofs / 3; row_bb++) {
269
270 t1(i) = w * row_diff_base_functions(i);
271
272 auto base_functions = col_data.getFTensor0N(gg, 0);
273 for (int col_bb = 0; col_bb != col_nb_dofs; col_bb++) {
274
275 FTensor::Tensor1<double *, 3> k(&locG(3 * row_bb + 0, col_bb),
276 &locG(3 * row_bb + 1, col_bb),
277 &locG(3 * row_bb + 2, col_bb));
278
279 k(i) += t1(i) * base_functions;
280 ++base_functions;
281 }
282 ++row_diff_base_functions;
283 }
284 for (; row_bb != row_nb_base_functions; row_bb++) {
285 ++row_diff_base_functions;
286 }
287 }
288
289 CHKERR MatSetValues(getFEMethod()->ksp_B, row_nb_dofs,
290 &*row_data.getIndices().begin(), col_nb_dofs,
291 &*col_data.getIndices().begin(), &*locG.data().begin(),
292 ADD_VALUES);
293
294 // ASSEMBLE THE TRANSPOSE
295 locG = trans(locG);
296 CHKERR MatSetValues(getFEMethod()->ksp_B, col_nb_dofs,
297 &*col_data.getIndices().begin(), row_nb_dofs,
298 &*row_data.getIndices().begin(), &*locG.data().begin(),
299 ADD_VALUES);
301 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
MoFEMErrorCode getBlockData(BlockData &data)

Member Data Documentation

◆ commonData

DataAtIntegrationPts& OpAssembleG::commonData

◆ dAta

BlockData& OpAssembleG::dAta

◆ locG

MatrixDouble OpAssembleG::locG

The documentation for this struct was generated from the following file: