v0.14.0
Public Member Functions | Public Attributes | List of all members
FractureMechanics::GriffithForceElement::OpLhs Struct Reference

#include <users_modules/fracture_mechanics/src/GriffithForceElement.hpp>

Inheritance diagram for FractureMechanics::GriffithForceElement::OpLhs:
[legend]
Collaboration diagram for FractureMechanics::GriffithForceElement::OpLhs:
[legend]

Public Member Functions

 OpLhs (int tag, BlockData &block_data, CommonData &common_data, const double beta_gc=0)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 
- Public Member Functions inherited from FractureMechanics::GriffithForceElement::AuxOp
 AuxOp (int tag, BlockData &block_data, CommonData &common_data)
 
MoFEMErrorCode setIndices (DataForcesAndSourcesCore::EntData &data)
 
MoFEMErrorCode setVariables (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, DataForcesAndSourcesCore::EntData &data)
 
MoFEMErrorCode setLambdaNodes (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, const std::string &lambda_field_name)
 
MoFEMErrorCode setLambdaIndices (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, const std::string &lambda_field_name)
 

Public Attributes

MatrixDouble mat
 
const double betaGC
 
- Public Attributes inherited from FractureMechanics::GriffithForceElement::AuxOp
int tAg
 
BlockDatablockData
 
CommonDatacommonData
 
ublas::vector< int > rowIndices
 
ublas::vector< int > rowIndicesLocal
 
VectorInt3 rowLambdaIndices
 
VectorInt3 rowLambdaIndicesLocal
 
VectorDouble3 lambdaAtNodes
 
VectorDouble activeVariables
 

Detailed Description

Definition at line 601 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ OpLhs()

FractureMechanics::GriffithForceElement::OpLhs::OpLhs ( int  tag,
BlockData block_data,
CommonData common_data,
const double  beta_gc = 0 
)
inline

Definition at line 604 of file GriffithForceElement.hpp.

607  "MESH_NODE_POSITIONS", "MESH_NODE_POSITIONS",
609  betaGC(beta_gc), AuxOp(tag, block_data, common_data) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode FractureMechanics::GriffithForceElement::OpLhs::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)
inline

Definition at line 614 of file GriffithForceElement.hpp.

617  {
619 
620  if (row_type != MBVERTEX) {
622  }
623 
631 
632  AuxFunctions<double> auxFun;
633  auxFun.currentCoords.resize(9, false);
634 
635  int nb_dofs = row_data.getFieldData().size();
636  for (int dd = 0; dd != nb_dofs; dd++)
637  auxFun.currentCoords[dd] = row_data.getFieldData()[dd];
638 
639  CHKERR setIndices(row_data);
640  CHKERR setVariables(this, row_data);
641  int nb_gauss_pts = row_data.getN().size1();
642  int nb_base_row = row_data.getFieldData().size() / 3;
643  int nb_base_col = col_data.getFieldData().size() / 3;
644  int row_nb_dofs = row_data.getIndices().size();
645 
646  mat.resize(9, 9, false);
647  mat.clear();
648  auto calculate_derivative = [&](FTensor::Tensor1<double *, 2> &t_diff) {
650  t_d_n(i, j) = FTensor::levi_civita(i, j, k) * auxFun.t_Tangent1_ksi(k) *
651  t_diff(N1);
652  t_d_n(i, j) -= FTensor::levi_civita(i, j, k) *
653  auxFun.t_Tangent2_eta(k) * t_diff(N0);
654 
655  return t_d_n;
656  };
657 
658  CHKERR auxFun.calculateAreaAndNormal(col_data.getDiffN());
659  auto t_normal = auxFun.t_Normal;
660  const double coefficient_1 = 0.5 * pow((t_normal(i) * t_normal(i)), -0.5);
661  const double coefficient_2 = 0.5 * pow((t_normal(i) * t_normal(i)), -1.5);
662 
663  // for homogeneous case, this is a bit ugly
664  auto &rho_v = commonData.densityRho;
665  if (rho_v.empty() || rho_v.size() != nb_base_row) {
666  rho_v.resize(nb_base_row, false);
667  }
668  auto rho = getFTensor0FromVec(rho_v);
669  for (int gg = 0; gg != nb_gauss_pts; gg++) {
670 
671  double val = getGaussPts()(2, gg) * 0.5;
672  auto t_row_diff = row_data.getFTensor1DiffN<2>(gg, 0);
673  for (int rrr = 0; rrr != nb_base_row; rrr++) {
674 
676  &mat(3 * rrr + 0, 0), &mat(3 * rrr + 0, 1), &mat(3 * rrr + 0, 2),
677  &mat(3 * rrr + 1, 0), &mat(3 * rrr + 1, 1), &mat(3 * rrr + 1, 2),
678  &mat(3 * rrr + 2, 0), &mat(3 * rrr + 2, 1), &mat(3 * rrr + 2, 2),
679  3);
680 
681  auto t_tan_row = calculate_derivative(t_row_diff);
682  auto t_col_diff = col_data.getFTensor1DiffN<2>(gg, 0);
683  const double coef = blockData.gc * val * pow(rho, betaGC);
684  for (int ccc = 0; ccc != nb_base_col; ccc++) {
685 
686  auto t_tan_col = calculate_derivative(t_col_diff);
687 
688  t_mat(i, j) -= coefficient_1 *
689  (t_tan_row(i, l) * t_tan_col(l, j) +
690  FTensor::levi_civita(i, j, k) * t_col_diff(N0) *
691  t_row_diff(N1) * t_normal(k) -
692  FTensor::levi_civita(i, j, k) * t_col_diff(N1) *
693  t_row_diff(N0) * t_normal(k));
694 
695  t_mat(i, j) += coefficient_2 * (t_tan_row(i, k) * t_normal(k)) *
696  (t_tan_col(l, j) * t_normal(l));
697 
698  t_mat(i, j) *= coef;
699  ++t_col_diff;
700  ++t_mat;
701  }
702  ++t_row_diff;
703  ++rho;
704  }
705  }
706 
707  int col_nb_dofs = col_data.getIndices().size();
708 
709  CHKERR MatSetValues(getFEMethod()->snes_B, row_nb_dofs,
710  &*rowIndices.data().begin(), col_nb_dofs,
711  &*col_data.getIndices().data().begin(),
712  &*mat.data().begin(), ADD_VALUES);
713 
715  }

Member Data Documentation

◆ betaGC

const double FractureMechanics::GriffithForceElement::OpLhs::betaGC

Definition at line 612 of file GriffithForceElement.hpp.

◆ mat

MatrixDouble FractureMechanics::GriffithForceElement::OpLhs::mat

Definition at line 611 of file GriffithForceElement.hpp.


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
FractureMechanics::GriffithForceElement::AuxOp::AuxOp
AuxOp(int tag, BlockData &block_data, CommonData &common_data)
Definition: GriffithForceElement.hpp:221
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
rho
double rho
Definition: plastic.cpp:140
FractureMechanics::GriffithForceElement::OpLhs::mat
MatrixDouble mat
Definition: GriffithForceElement.hpp:611
FractureMechanics::GriffithForceElement::OpLhs::betaGC
const double betaGC
Definition: GriffithForceElement.hpp:612
FractureMechanics::GriffithForceElement::AuxOp::setIndices
MoFEMErrorCode setIndices(DataForcesAndSourcesCore::EntData &data)
Definition: GriffithForceElement.hpp:232
FTensor::levi_civita
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
Definition: Levi_Civita.hpp:617
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
FractureMechanics::GriffithForceElement::BlockData::gc
double gc
Griffith energy.
Definition: GriffithForceElement.hpp:75
FTensor::Tensor2< double, 3, 3 >
FractureMechanics::GriffithForceElement::AuxOp::setVariables
MoFEMErrorCode setVariables(FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, DataForcesAndSourcesCore::EntData &data)
Definition: GriffithForceElement.hpp:261
FTensor::Number< 0 >
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
FractureMechanics::GriffithForceElement::AuxOp::blockData
BlockData & blockData
Definition: GriffithForceElement.hpp:218
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
FractureMechanics::GriffithForceElement::CommonData::densityRho
VectorDouble densityRho
gc * rho^beta
Definition: GriffithForceElement.hpp:69
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
FractureMechanics::GriffithForceElement::AuxOp::rowIndices
ublas::vector< int > rowIndices
Definition: GriffithForceElement.hpp:222
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:429
FractureMechanics::GriffithForceElement::AuxOp::commonData
CommonData & commonData
Definition: GriffithForceElement.hpp:219
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21