v0.14.0
Public Member Functions | Public Attributes | List of all members
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 *.

Examples
elasticity_mixed_formulation.cpp, and ElasticityMixedFormulation.hpp.

Definition at line 215 of file ElasticityMixedFormulation.hpp.

Constructor & Destructor Documentation

◆ OpAssembleG()

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

Definition at line 222 of file ElasticityMixedFormulation.hpp.

224  "U", "P", UserDataOperator::OPROWCOL),
225  commonData(common_data), dAta(data) {
226  sYmm = false;
227  }

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 229 of file ElasticityMixedFormulation.hpp.

232  {
233 
235 
236  const int row_nb_dofs = row_data.getIndices().size();
237  if (!row_nb_dofs)
239  const int col_nb_dofs = col_data.getIndices().size();
240  if (!col_nb_dofs)
242 
243  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
244  dAta.tEts.end()) {
246  }
248 
249  // Set size can clear local tangent matrix
250  locG.resize(row_nb_dofs, col_nb_dofs, false);
251  locG.clear();
252  const int row_nb_gauss_pts = row_data.getN().size1();
253  if (!row_nb_gauss_pts)
255  const int row_nb_base_functions = row_data.getN().size2();
256  auto row_diff_base_functions = row_data.getFTensor1DiffN<3>();
257 
261 
262  // INTEGRATION
263  for (int gg = 0; gg != row_nb_gauss_pts; gg++) {
264 
265  // Get volume and integration weight
266  double w = getVolume() * getGaussPts()(3, gg);
267 
268  int row_bb = 0;
269  for (; row_bb != row_nb_dofs / 3; row_bb++) {
270 
271  t1(i) = w * row_diff_base_functions(i);
272 
273  auto base_functions = col_data.getFTensor0N(gg, 0);
274  for (int col_bb = 0; col_bb != col_nb_dofs; col_bb++) {
275 
276  FTensor::Tensor1<double *, 3> k(&locG(3 * row_bb + 0, col_bb),
277  &locG(3 * row_bb + 1, col_bb),
278  &locG(3 * row_bb + 2, col_bb));
279 
280  k(i) += t1(i) * base_functions;
281  ++base_functions;
282  }
283  ++row_diff_base_functions;
284  }
285  for (; row_bb != row_nb_base_functions; row_bb++) {
286  ++row_diff_base_functions;
287  }
288  }
289 
290  CHKERR MatSetValues(getFEMethod()->ksp_B, row_nb_dofs,
291  &*row_data.getIndices().begin(), col_nb_dofs,
292  &*col_data.getIndices().begin(), &*locG.data().begin(),
293  ADD_VALUES);
294 
295  // ASSEMBLE THE TRANSPOSE
296  locG = trans(locG);
297  CHKERR MatSetValues(getFEMethod()->ksp_B, col_nb_dofs,
298  &*col_data.getIndices().begin(), row_nb_dofs,
299  &*row_data.getIndices().begin(), &*locG.data().begin(),
300  ADD_VALUES);
302  }

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:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1644
FTensor::Tensor1< double, 3 >
DataAtIntegrationPts::getBlockData
MoFEMErrorCode getBlockData(BlockData &data)
Definition: ElasticityMixedFormulation.hpp:52
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
OpAssembleG::locG
MatrixDouble locG
Definition: ElasticityMixedFormulation.hpp:219
OpAssembleG::commonData
DataAtIntegrationPts & commonData
Definition: ElasticityMixedFormulation.hpp:218
OpAssembleG::dAta
BlockData & dAta
Definition: ElasticityMixedFormulation.hpp:220
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
BlockData::tEts
Range tEts
Definition: ElasticityMixedFormulation.hpp:17
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
sdf_wavy_2d.w
int w
Definition: sdf_wavy_2d.py:6
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
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