v0.15.4
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
EshelbianPlasticity::OpCohesive_dJ_dkappa Struct Reference
Inheritance diagram for EshelbianPlasticity::OpCohesive_dJ_dkappa:
[legend]
Collaboration diagram for EshelbianPlasticity::OpCohesive_dJ_dkappa:
[legend]

Public Types

using OP = OpCohesiveRhs
 
- Public Types inherited from EshelbianPlasticity::OpCohesiveRhs
using OP = OpBrokenBaseCohesive
 
- Public Types inherited from EshelbianPlasticity::OpBrokenBaseCohesive
using OP = FormsIntegrators< FaceElementForcesAndSourcesCore::UserDataOperator >::Assembly< A >::OpBrokenBase
 
- Public Types inherited from MoFEM::OpBrokenBaseImpl< OpBase >
using OP = OpBase
 
- Public Types inherited from MoFEM::OpBaseImpl< A, EleOp >
using OpType = typename EleOp::OpType
 
using EntData = EntitiesFieldData::EntData
 
using MatSetValuesHook = boost::function< MoFEMErrorCode(ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)>
 

Public Member Functions

MoFEMErrorCode iNtegrate (EntData &data)
 Class dedicated to integrate operator.
 
MoFEMErrorCode aSsemble (EntData &data)
 Assemble local vector into global vector.
 
 OpCohesiveRhs (boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< double > gc_ptr, boost::shared_ptr< VectorDouble > kappa_ptr, boost::shared_ptr< VectorDouble > kappa_delta_ptr, boost::shared_ptr< std::array< MatrixDouble, 2 > > lambda_ptr=nullptr, Tag dissipation_tags=Tag(), Tag grad_dissipation_tags=Tag(), SmartPetscObj< Vec > vec_dJ_dx=SmartPetscObj< Vec >(), boost::shared_ptr< Range > ents_ptr=nullptr)
 
- Public Member Functions inherited from EshelbianPlasticity::OpCohesiveRhs
 OpCohesiveRhs (boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< double > gc_ptr, boost::shared_ptr< VectorDouble > kappa_ptr, boost::shared_ptr< VectorDouble > kappa_delta_ptr, boost::shared_ptr< std::array< MatrixDouble, 2 > > lambda_ptr=nullptr, Tag dissipation_tags=Tag(), Tag grad_dissipation_tags=Tag(), SmartPetscObj< Vec > vec_dJ_dx=SmartPetscObj< Vec >(), boost::shared_ptr< Range > ents_ptr=nullptr)
 
MoFEMErrorCode iNtegrate (EntData &data)
 Class dedicated to integrate operator.
 
- Public Member Functions inherited from EshelbianPlasticity::OpBrokenBaseCohesive
 OpBrokenBaseCohesive (boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_flux_data_ptr, boost::shared_ptr< Range > ents_ptr=nullptr)
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 
- Public Member Functions inherited from MoFEM::OpBrokenBaseImpl< OpBase >
 OpBrokenBaseImpl (boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< Range > ents_ptr=nullptr)
 
 OpBrokenBaseImpl (const std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, const bool assmb_transpose, const bool only_transpose, boost::shared_ptr< Range > ents_ptr=nullptr)
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 
- Public Member Functions inherited from MoFEM::OpBaseImpl< A, EleOp >
 OpBaseImpl (const std::string row_field_name, const std::string col_field_name, const OpType type, boost::shared_ptr< Range > ents_ptr=nullptr)
 Constructor for base operator implementation.
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Do calculations for the left hand side.
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntData &row_data)
 Do calculations for the right hand side.
 

Additional Inherited Members

- Public Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
TimeFun timeScalingFun
 assumes that time variable is set
 
FEFun feScalingFun
 set by fe entity handle
 
boost::shared_ptr< RangeentsPtr
 Entities on which element is run.
 
- Static Public Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
static MatSetValuesHook matSetValuesHook
 
- Protected Member Functions inherited from MoFEM::OpBaseImpl< A, EleOp >
template<int DIM>
FTensor::Tensor1< FTensor::PackPtr< double *, DIM >, DIM > getNf ()
 Get local vector tensor for assembly.
 
template<int DIM>
FTensor::Tensor2< FTensor::PackPtr< double *, DIM >, DIM, DIM > getLocMat (const int rr)
 Get local matrix tensor for assembly.
 
virtual MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Integrate grad-grad operator.
 
virtual MoFEMErrorCode aSsemble (EntData &row_data, EntData &col_data, const bool trans)
 Assemble local matrix into global matrix.
 
virtual size_t getNbOfBaseFunctions (EntitiesFieldData::EntData &data)
 Get number of base functions.
 
- Protected Attributes inherited from EshelbianPlasticity::OpCohesiveRhs
boost::shared_ptr< MatrixDouble > uGammaPtr
 
boost::shared_ptr< doublegcPtr
 
boost::shared_ptr< VectorDouble > kappaPtr
 
boost::shared_ptr< VectorDouble > kappaDeltaPtr
 
boost::shared_ptr< std::array< MatrixDouble, 2 > > lambdaPtr
 
boost::shared_ptr< doubletotalDissipation
 
boost::shared_ptr< doubletatalDissipationGrad
 
SmartPetscObj< Vec > vec_dJdu
 
Tag dissipationTag
 
Tag gradDissipationTag
 
- Protected Attributes inherited from EshelbianPlasticity::OpBrokenBaseCohesive
boost::shared_ptr< MatrixDouble > fluxMatPtr
 
int faceSense = 0
 
- Protected Attributes inherited from MoFEM::OpBrokenBaseImpl< OpBase >
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
 
- Protected Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
int nbRows
 number of dofs on rows
 
int nbCols
 number if dof on column
 
int nbIntegrationPts
 number of integration points
 
int nbRowBaseFunctions
 number or row base functions
 
int rowSide
 row side number
 
int colSide
 column side number
 
EntityType rowType
 row type
 
EntityType colType
 column type
 
bool assembleTranspose
 
bool onlyTranspose
 
MatrixDouble locMat
 local entity block matrix
 
MatrixDouble locMatTranspose
 local entity block matrix
 
VectorDouble locF
 local entity vector
 

Detailed Description

Definition at line 565 of file EshelbianCohesive.cpp.

Member Typedef Documentation

◆ OP

Definition at line 567 of file EshelbianCohesive.cpp.

Member Function Documentation

◆ aSsemble()

MoFEMErrorCode EshelbianPlasticity::OpCohesive_dJ_dkappa::aSsemble ( EntData data)
inlinevirtual

Assemble local vector into global vector.

Parameters
dataEntity data
Returns
Error code

Reimplemented from MoFEM::OpBaseImpl< A, EleOp >.

Definition at line 666 of file EshelbianCohesive.cpp.

666 {
669 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()

◆ iNtegrate()

MoFEMErrorCode EshelbianPlasticity::OpCohesive_dJ_dkappa::iNtegrate ( EntData data)
inlinevirtual

Class dedicated to integrate operator.

Parameters
dataentity data on element row
Returns
error code

Reimplemented from MoFEM::OpBaseImpl< A, EleOp >.

Definition at line 570 of file EshelbianCohesive.cpp.

570 {
572
573 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
574
575 FTENSOR_INDEX(3, i);
576 FTENSOR_INDEX(3, J);
577 FTENSOR_INDEX(3, K);
578
579 int nb_integration_pts = OP::getGaussPts().size2();
580
581 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
582 auto t_lambda = getFTensor2FromMat<3, 3>(lambdaPtr->at(get_sense_index()));
583 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
584 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
585 auto t_face_normal = getFTensor1NormalsAtGaussPts();
586 auto t_w = OP::getFTensor0IntegrationWeight();
587
588 auto next = [&]() {
589 ++t_P;
590 ++t_face_normal;
591 ++t_kappa;
592 ++t_delta_kappa;
593 ++t_lambda;
594 ++t_w;
595 };
596
597 double face_dissipation = 0.0;
598 double face_grad_dissipation = 0.0;
599 auto face_handle = OP::getFEEntityHandle();
600 CHKERR getMoab().tag_get_data(dissipationTag, &face_handle, 1,
601 &face_dissipation);
602 CHKERR getMoab().tag_get_data(gradDissipationTag, &face_handle, 1,
603 &face_grad_dissipation);
604
605 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
607 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
608 FTensor::Tensor1<double, 3> t_normalized_normal;
609 t_normalized_normal(J) = t_normal(J);
610 t_normalized_normal.normalize();
612 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
613
614 auto dJ_dkappa = [](auto &t_delta_kappa, auto &t_traction,
615 auto &t_normalized_normal, auto &t_kappa, auto gc) {
616 double kappa = static_cast<double>(t_kappa);
617 double delta_kappa = static_cast<double>(t_delta_kappa);
619 t_traction, t_normalized_normal, OpGetParameters::beta);
621 delta_kappa, teff, kappa, gc);
622 double m_grad =
624 delta_kappa, teff, kappa, gc);
625 return boost::make_tuple(m, m_grad);
626 };
627
628 auto dr_kappa = [](auto &t_delta_kappa, auto &t_traction,
629 auto &t_normalized_normal, auto &t_kappa, auto gc) {
630 double kappa = static_cast<double>(t_kappa);
631 double delta_kappa = static_cast<double>(t_delta_kappa);
632 double kappa_plus_delta = kappa + delta_kappa;
633 double diff_alpha =
634 GriffithCohesiveLaw::getDiffAlpha(kappa_plus_delta, gc);
636 t_traction, t_normalized_normal, diff_alpha,
638 FTENSOR_INDEX(3, i);
639 FTensor::Tensor1<double, 3> t_gap_double;
640 t_gap_double(i) = -t_gap(i);
641 return t_gap_double;
642 };
643
644 auto [J, dJ] = dJ_dkappa(t_delta_kappa, t_traction, t_normalized_normal,
645 t_kappa, *gcPtr);
646 face_dissipation += t_w * J * t_normal.l2();
647 face_grad_dissipation += t_w * dJ * t_normal.l2();
648
649 auto t_dr = dr_kappa(t_delta_kappa, t_traction, t_normalized_normal,
650 t_kappa, *gcPtr);
652 t_l(i) = t_lambda(i, K) * t_normal(K);
653 face_grad_dissipation -= t_w * t_l(i) * t_dr(i);
654
655 next();
656 }
657
658 CHKERR getMoab().tag_set_data(dissipationTag, &face_handle, 1,
659 &face_dissipation);
660 CHKERR getMoab().tag_set_data(gradDissipationTag, &face_handle, 1,
661 &face_grad_dissipation);
662
664 };
#define FTENSOR_INDEX(DIM, I)
Tensor1< T, Tensor_Dim > normalize()
#define CHKERR
Inline error check.
double kappa
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'm', 3 > m
static auto calculateDissipationSurplus(const T &delta_kappa, const T t_eff, const T &kappa, double gc)
static auto calculateGap(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double alpha, double beta)
static auto calculateEffectiveTraction(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double beta)
static auto calculateDissipationSurplusDiffKappa(const T &delta_kappa, const T &t_eff, const T &kappa, double gc)
static T getDiffAlpha(const T &kappa, double gc)
boost::shared_ptr< MatrixDouble > fluxMatPtr
boost::shared_ptr< double > gcPtr
boost::shared_ptr< VectorDouble > kappaDeltaPtr
boost::shared_ptr< VectorDouble > kappaPtr
boost::shared_ptr< std::array< MatrixDouble, 2 > > lambdaPtr

◆ OpCohesiveRhs()

EshelbianPlasticity::OpCohesiveRhs::OpCohesiveRhs ( boost::shared_ptr< std::vector< BrokenBaseSideData > >  broken_base_side_data,
boost::shared_ptr< double gc_ptr,
boost::shared_ptr< VectorDouble >  kappa_ptr,
boost::shared_ptr< VectorDouble >  kappa_delta_ptr,
boost::shared_ptr< std::array< MatrixDouble, 2 > >  lambda_ptr = nullptr,
Tag  dissipation_tags = Tag(),
Tag  grad_dissipation_tags = Tag(),
SmartPetscObj< Vec >  vec_dJ_dx = SmartPetscObj<Vec>(),
boost::shared_ptr< Range ents_ptr = nullptr 
)
inline

Definition at line 366 of file EshelbianCohesive.cpp.

375 : OpBrokenBaseCohesive(broken_base_side_data, ents_ptr), gcPtr(gc_ptr),
376 kappaPtr(kappa_ptr), kappaDeltaPtr(kappa_delta_ptr),
377 lambdaPtr(lambda_ptr), dissipationTag(dissipation_tags),
378 gradDissipationTag(grad_dissipation_tags), vec_dJdu(vec_dJ_dx) {}
OpBrokenBaseCohesive(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_flux_data_ptr, boost::shared_ptr< Range > ents_ptr=nullptr)

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