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

Assemble P *. More...

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

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

Public Member Functions

 OpAssembleP (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 locP
 
MatrixDouble translocP
 
BlockDatadAta
 

Detailed Description

Assemble P *.

Examples
elasticity_mixed_formulation.cpp, and ElasticityMixedFormulation.hpp.

Definition at line 114 of file ElasticityMixedFormulation.hpp.

Constructor & Destructor Documentation

◆ OpAssembleP()

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

Definition at line 122 of file ElasticityMixedFormulation.hpp.

124  "P", "P", UserDataOperator::OPROWCOL),
125  commonData(common_data), dAta(data) {
126  sYmm = true;
127  }

Member Function Documentation

◆ doWork()

PetscErrorCode OpAssembleP::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 129 of file ElasticityMixedFormulation.hpp.

132  {
133 
135  const int row_nb_dofs = row_data.getIndices().size();
136  if (!row_nb_dofs)
138  const int col_nb_dofs = col_data.getIndices().size();
139  if (!col_nb_dofs)
141 
142  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
143  dAta.tEts.end()) {
145  }
147  // Set size can clear local tangent matrix
148  locP.resize(row_nb_dofs, col_nb_dofs, false);
149  locP.clear();
150 
151  const int row_nb_gauss_pts = row_data.getN().size1();
152  if (!row_nb_gauss_pts)
154  const int row_nb_base_functions = row_data.getN().size2();
155  auto row_base_functions = row_data.getFTensor0N();
156 
157  // get data
160  const double lambda = commonData.lAmbda;
161 
162  double coefficient = commonData.pOisson == 0.5 ? 0. : 1 / lambda;
163 
164  // integration
165  if (coefficient != 0.) {
166  for (int gg = 0; gg != row_nb_gauss_pts; gg++) {
167 
168  // Get volume and integration weight
169  double w = getVolume() * getGaussPts()(3, gg);
170 
171  // INTEGRATION
172  int row_bb = 0;
173  for (; row_bb != row_nb_dofs; row_bb++) {
174  auto col_base_functions = col_data.getFTensor0N(gg, 0);
175  for (int col_bb = 0; col_bb != col_nb_dofs; col_bb++) {
176 
177  locP(row_bb, col_bb) -=
178  w * row_base_functions * col_base_functions * coefficient;
179 
180  ++col_base_functions;
181  }
182  ++row_base_functions;
183  }
184  for (; row_bb != row_nb_base_functions; row_bb++) {
185  ++row_base_functions;
186  }
187  }
188  }
189 
190 
192  getFEMethod()->ksp_B, row_nb_dofs, &*row_data.getIndices().begin(),
193  col_nb_dofs, &*col_data.getIndices().begin(), &locP(0, 0), ADD_VALUES);
194 
195  // is symmetric
196  if (row_side != col_side || row_type != col_type) {
197  translocP.resize(col_nb_dofs, row_nb_dofs, false);
198  noalias(translocP) = trans(locP);
199  CHKERR MatSetValues(getFEMethod()->ksp_B, col_nb_dofs,
200  &*col_data.getIndices().begin(), row_nb_dofs,
201  &*row_data.getIndices().begin(), &translocP(0, 0),
202  ADD_VALUES);
203  }
204 
206  }

Member Data Documentation

◆ commonData

DataAtIntegrationPts& OpAssembleP::commonData

◆ dAta

BlockData& OpAssembleP::dAta

◆ locP

MatrixDouble OpAssembleP::locP

◆ translocP

MatrixDouble OpAssembleP::translocP

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
PlasticOps::w
double w(double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y)
Definition: PlasticOpsGeneric.hpp:276
DataAtIntegrationPts::pOisson
double pOisson
Definition: ElasticityMixedFormulation.hpp:26
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
DataAtIntegrationPts::lAmbda
double lAmbda
Definition: ElasticityMixedFormulation.hpp:28
OpAssembleP::translocP
MatrixDouble translocP
Definition: ElasticityMixedFormulation.hpp:119
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
OpAssembleP::commonData
DataAtIntegrationPts & commonData
Definition: ElasticityMixedFormulation.hpp:117
BlockData::tEts
Range tEts
Definition: ElasticityMixedFormulation.hpp:17
OpAssembleP::dAta
BlockData & dAta
Definition: ElasticityMixedFormulation.hpp:120
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
OpAssembleP::locP
MatrixDouble locP
Definition: ElasticityMixedFormulation.hpp:118
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