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

Assemble K *. More...

#include <users_modules/tutorials/cor-7/src/ElasticityMixedFormulation.hpp>

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

Public Member Functions

 OpAssembleK (DataAtIntegrationPts &common_data, BlockData &data, bool symm=true)
 
PetscErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 

Public Attributes

MatrixDouble locK
 
MatrixDouble translocK
 
FTensor::Tensor2< double, 3, 3 > diffDiff
 
DataAtIntegrationPtscommonData
 
BlockDatadAta
 

Detailed Description

Assemble K *.

Examples
elasticity_mixed_formulation.cpp, and ElasticityMixedFormulation.hpp.

Definition at line 311 of file ElasticityMixedFormulation.hpp.

Constructor & Destructor Documentation

◆ OpAssembleK()

OpAssembleK::OpAssembleK ( DataAtIntegrationPts common_data,
BlockData data,
bool  symm = true 
)
inline
Examples
ElasticityMixedFormulation.hpp.

Definition at line 321 of file ElasticityMixedFormulation.hpp.

324  symm),
325  commonData(common_data), dAta(data) {}

Member Function Documentation

◆ doWork()

PetscErrorCode OpAssembleK::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 327 of file ElasticityMixedFormulation.hpp.

330  {
332 
333  auto get_tensor2 = [](MatrixDouble &m, const int r, const int c) {
335  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2),
336  &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2),
337  &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2));
338  };
339 
340  const int row_nb_dofs = row_data.getIndices().size();
341  if (!row_nb_dofs)
343  const int col_nb_dofs = col_data.getIndices().size();
344  if (!col_nb_dofs)
346  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
347  dAta.tEts.end()) {
349  }
351 
352  const bool diagonal_block =
353  (row_type == col_type) && (row_side == col_side);
354  // get number of integration points
355  // Set size can clear local tangent matrix
356  locK.resize(row_nb_dofs, col_nb_dofs, false);
357  locK.clear();
358 
359  const int row_nb_gauss_pts = row_data.getN().size1();
360  const int row_nb_base_functions = row_data.getN().size2();
361 
362  auto row_diff_base_functions = row_data.getFTensor1DiffN<3>();
363 
368 
369  // integrate local matrix for entity block
370  for (int gg = 0; gg != row_nb_gauss_pts; gg++) {
371 
372  // Get volume and integration weight
373  double w = getVolume() * getGaussPts()(3, gg);
374 
375  int row_bb = 0;
376  for (; row_bb != row_nb_dofs / 3; row_bb++) {
377 
378  auto col_diff_base_functions = col_data.getFTensor1DiffN<3>(gg, 0);
379  const int final_bb = diagonal_block ? row_bb + 1 : col_nb_dofs / 3;
380  int col_bb = 0;
381  for (; col_bb != final_bb; col_bb++) {
382 
383  auto t_assemble = get_tensor2(locK, 3 * row_bb, 3 * col_bb);
384 
385  diffDiff(j, l) =
386  w * row_diff_base_functions(j) * col_diff_base_functions(l);
387 
388  t_assemble(i, k) += diffDiff(j, l) * commonData.tD(i, j, k, l);
389  // Next base function for column
390  ++col_diff_base_functions;
391  }
392 
393  ++row_diff_base_functions;
394  }
395  for (; row_bb != row_nb_base_functions; row_bb++) {
396  ++row_diff_base_functions;
397  }
398  }
399 
400  if (diagonal_block) {
401  for (int row_bb = 0; row_bb != row_nb_dofs / 3; row_bb++) {
402  int col_bb = 0;
403  for (; col_bb != row_bb + 1; col_bb++) {
404  auto t_assemble = get_tensor2(locK, 3 * row_bb, 3 * col_bb);
405  auto t_off_side = get_tensor2(locK, 3 * col_bb, 3 * row_bb);
406  t_off_side(i, k) = t_assemble(k, i);
407  }
408  }
409  }
410 
411  const int *row_ind = &*row_data.getIndices().begin();
412  const int *col_ind = &*col_data.getIndices().begin();
413  Mat B = getFEMethod()->ksp_B != PETSC_NULL ? getFEMethod()->ksp_B
414  : getFEMethod()->ksp_B;
415  CHKERR MatSetValues(B, row_nb_dofs, row_ind, col_nb_dofs, col_ind,
416  &locK(0, 0), ADD_VALUES);
417 
418  if (row_type != col_type || row_side != col_side) {
419  translocK.resize(col_nb_dofs, row_nb_dofs, false);
420  noalias(translocK) = trans(locK);
421  CHKERR MatSetValues(B, col_nb_dofs, col_ind, row_nb_dofs, row_ind,
422  &translocK(0, 0), ADD_VALUES);
423  }
424 
426  }

Member Data Documentation

◆ commonData

DataAtIntegrationPts& OpAssembleK::commonData

◆ dAta

BlockData& OpAssembleK::dAta

◆ diffDiff

FTensor::Tensor2<double, 3, 3> OpAssembleK::diffDiff

◆ locK

MatrixDouble OpAssembleK::locK

◆ translocK

MatrixDouble OpAssembleK::translocK

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::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:1631
OpAssembleK::translocK
MatrixDouble translocK
Definition: ElasticityMixedFormulation.hpp:315
PlasticOps::w
double w(double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y)
Definition: PlasticOpsGeneric.hpp:276
OpAssembleK::dAta
BlockData & dAta
Definition: ElasticityMixedFormulation.hpp:319
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
DataAtIntegrationPts::getBlockData
MoFEMErrorCode getBlockData(BlockData &data)
Definition: ElasticityMixedFormulation.hpp:52
sdf.r
int r
Definition: sdf.py:8
FTensor::Tensor2< double *, 3, 3 >
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
OpAssembleK::diffDiff
FTensor::Tensor2< double, 3, 3 > diffDiff
Definition: ElasticityMixedFormulation.hpp:316
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
OpAssembleK::locK
MatrixDouble locK
Definition: ElasticityMixedFormulation.hpp:314
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
DataAtIntegrationPts::tD
FTensor::Ddg< double, 3, 3 > tD
Definition: ElasticityMixedFormulation.hpp:24
OpAssembleK::commonData
DataAtIntegrationPts & commonData
Definition: ElasticityMixedFormulation.hpp:318
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
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:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21