v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
OpCalculateMassMatrix Struct Reference

Assemble mass matrix for elastic element TODO: CHANGE FORMULA *. More...

#include <users_modules/basic_finite_elements/src/HookeElement.hpp>

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

Public Member Functions

 OpCalculateMassMatrix (const std::string row_field, const std::string col_field, BlockData &data, MassBlockData &mass_data, boost::shared_ptr< DataAtIntegrationPts > &common_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
 
BlockDatadAta
 
MassBlockDatamassData
 
boost::shared_ptr< DataAtIntegrationPtscommonData
 

Detailed Description

Assemble mass matrix for elastic element TODO: CHANGE FORMULA *.

*

\[ {\bf{M}} = \int\limits_\Omega \]

Examples
HookeElement.hpp.

Definition at line 212 of file HookeElement.hpp.

Constructor & Destructor Documentation

◆ OpCalculateMassMatrix()

OpCalculateMassMatrix::OpCalculateMassMatrix ( const std::string  row_field,
const std::string  col_field,
BlockData data,
MassBlockData mass_data,
boost::shared_ptr< DataAtIntegrationPts > &  common_data,
bool  symm = true 
)
inline

Definition at line 222 of file HookeElement.hpp.

227 : VolumeElementForcesAndSourcesCore::UserDataOperator(
228 row_field, col_field, OPROWCOL, symm),
229 commonData(common_data), dAta(data), massData(mass_data) {}
MassBlockData & massData
boost::shared_ptr< DataAtIntegrationPts > commonData

Member Function Documentation

◆ doWork()

PetscErrorCode OpCalculateMassMatrix::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
inline

Definition at line 231 of file HookeElement.hpp.

234 {
236
237 auto get_tensor2 = [](MatrixDouble &m, const int r, const int c) {
239 &m(3 * r + 0, 3 * c + 0), &m(3 * r + 0, 3 * c + 1),
240 &m(3 * r + 0, 3 * c + 2), &m(3 * r + 1, 3 * c + 0),
241 &m(3 * r + 1, 3 * c + 1), &m(3 * r + 1, 3 * c + 2),
242 &m(3 * r + 2, 3 * c + 0), &m(3 * r + 2, 3 * c + 1),
243 &m(3 * r + 2, 3 * c + 2));
244 };
245
246 const int row_nb_dofs = row_data.getIndices().size();
247 if (!row_nb_dofs)
249 const int col_nb_dofs = col_data.getIndices().size();
250 if (!col_nb_dofs)
252 if (dAta.tEts.find(getFEEntityHandle()) == dAta.tEts.end()) {
254 }
255 if (massData.tEts.find(getFEEntityHandle()) == massData.tEts.end()) {
257 }
258
259 const bool diagonal_block =
260 (row_type == col_type) && (row_side == col_side);
261 // get number of integration points
262 // Set size can clear local tangent matrix
263 locK.resize(row_nb_dofs, col_nb_dofs, false);
264 locK.clear();
265
266 const int row_nb_gauss_pts = row_data.getN().size1();
267 const int row_nb_base_functions = row_data.getN().size2();
268
273
274 double density = massData.rho0;
275
276 // get integration weights
277 auto t_w = getFTensor0IntegrationWeight();
278
279 // integrate local matrix for entity block
280 for (int gg = 0; gg != row_nb_gauss_pts; gg++) {
281
282 auto t_row_base_func = row_data.getFTensor0N(gg, 0);
283
284 // Get volume and integration weight
285 double w = getVolume() * t_w;
286
287 for (int row_bb = 0; row_bb != row_nb_dofs / 3; row_bb++) {
288 auto t_col_base_func = col_data.getFTensor0N(gg, 0);
289 for (int col_bb = 0; col_bb != col_nb_dofs / 3; col_bb++) {
290 auto t_assemble = get_tensor2(locK, row_bb, col_bb);
291 t_assemble(i, j) += density * t_row_base_func * t_col_base_func * w;
292 // Next base function for column
293 ++t_col_base_func;
294 }
295 // Next base function for row
296 ++t_row_base_func;
297 }
298 // Next integration point for getting weight
299 ++t_w;
300 }
301
302 CHKERR MatSetValues(getKSPB(), row_data, col_data, &locK(0, 0),
303 ADD_VALUES);
304
305 // is symmetric
306 if (row_type != col_type || row_side != col_side) {
307 translocK.resize(col_nb_dofs, row_nb_dofs, false);
308 noalias(translocK) = trans(locK);
309
310 CHKERR MatSetValues(getKSPB(), col_data, row_data, &translocK(0, 0),
311 ADD_VALUES);
312 }
313
315 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'm', SPACE_DIM > m
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'l', 3 > l
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.
double w(const double g, const double t)
int r
Definition: sdf.py:8
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.

Member Data Documentation

◆ commonData

boost::shared_ptr<DataAtIntegrationPts> OpCalculateMassMatrix::commonData

Definition at line 220 of file HookeElement.hpp.

◆ dAta

BlockData& OpCalculateMassMatrix::dAta

Definition at line 217 of file HookeElement.hpp.

◆ locK

MatrixDouble OpCalculateMassMatrix::locK

Definition at line 215 of file HookeElement.hpp.

◆ massData

MassBlockData& OpCalculateMassMatrix::massData

Definition at line 218 of file HookeElement.hpp.

◆ translocK

MatrixDouble OpCalculateMassMatrix::translocK

Definition at line 216 of file HookeElement.hpp.


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