v0.14.0
Public Member Functions | Private Attributes | List of all members
EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP Struct Reference

#include <users_modules/eshelbian_plasticit/src/EshelbianContact.hpp>

Inheritance diagram for EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP:
[legend]
Collaboration diagram for EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP:
[legend]

Public Member Functions

 OpConstrainBoundaryL2Lhs_dP (const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< ContactOps::CommonData > common_data_ptr, boost::shared_ptr< ContactTree > contact_tree_ptr, boost::shared_ptr< std::map< int, Range >> sdf_map_range_ptr=nullptr)
 
MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 

Private Attributes

boost::shared_ptr< ContactOps::CommonDatacommonDataPtr
 
boost::shared_ptr< ContactTreecontactTreePtr
 
boost::shared_ptr< std::map< int, Range > > sdfMapRangePtr
 

Detailed Description

Examples
EshelbianPlasticity.cpp.

Definition at line 160 of file EshelbianContact.hpp.

Constructor & Destructor Documentation

◆ OpConstrainBoundaryL2Lhs_dP()

EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::OpConstrainBoundaryL2Lhs_dP ( const std::string  row_field_name,
const std::string  col_field_name,
boost::shared_ptr< ContactOps::CommonData common_data_ptr,
boost::shared_ptr< ContactTree contact_tree_ptr,
boost::shared_ptr< std::map< int, Range >>  sdf_map_range_ptr = nullptr 
)

Definition at line 638 of file EshelbianContact.cpp.

643  : ContactOps::AssemblyBoundaryEleOp(row_field_name, col_field_name,
645  commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
646  sdfMapRangePtr(sdf_map_range_ptr) {
647 
648  sYmm = false;
649 }

Member Function Documentation

◆ iNtegrate()

MoFEMErrorCode EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::iNtegrate ( EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)

Definition at line 652 of file EshelbianContact.cpp.

653  {
655 
656  using namespace ContactOps;
657 
661 
662  auto nb_rows = row_data.getIndices().size();
663  auto nb_cols = col_data.getIndices().size();
664 
665  auto &locMat = AssemblyBoundaryEleOp::locMat;
666  locMat.resize(nb_rows, nb_cols, false);
667  locMat.clear();
668 
669  if (nb_cols && nb_rows) {
670 
671  const size_t nb_gauss_pts = getGaussPts().size2();
672 
673  auto t_w = getFTensor0IntegrationWeight();
674  auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
675  auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
676  auto t_coords = getFTensor1CoordsAtGaussPts();
677  auto t_material_normal = getFTensor1NormalsAtGaussPts();
678 
679  // placeholder to pass boundary block id to python
680  auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
681  getSdf(this, commonDataPtr->contactDisp,
682  checkSdf(getFEEntityHandle(), *sdfMapRangePtr), false);
683 
684  auto t_sdf_v = getFTensor0FromVec(v_sdf);
685  auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
686 
687  auto next = [&]() {
688  ++t_w;
689  ++t_disp;
690  ++t_traction;
691  ++t_coords;
692  ++t_material_normal;
693  ++t_sdf_v;
694  ++t_grad_sdf_v;
695  };
696 
697  auto face_data_vec_ptr =
698  contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
699  auto face_gauss_pts_it = face_data_vec_ptr->begin();
700 
701  auto t_row_base = row_data.getFTensor0N();
702  auto nb_face_functions = row_data.getN().size2();
703 
704  constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
705 
706  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
707 
708  auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
709  face_data_vec_ptr);
710 
711  auto check_face_contact = [&]() {
712  if (face_data_ptr) {
713  return true;
714  }
715  return false;
716  };
717 
719 
720  auto tn = t_traction(i) * t_grad_sdf_v(i);
721  auto c = ContactOps::constrain(t_sdf_v, tn);
722  if (!c && check_face_contact()) {
723  FTensor::Tensor1<double, 3> t_spatial_coords;
724  t_spatial_coords(i) = t_coords(i) + t_disp(i);
725  constexpr double eps = std::numeric_limits<float>::epsilon();
726  for (auto ii = 0; ii != 3; ++ii) {
727  FTensor::Tensor1<std::complex<double>, 3> t_traction_cx{
728  t_traction(0), t_traction(1), t_traction(2)};
729  t_traction_cx(ii) += eps * 1i;
730  auto t_rhs_tmp =
731  multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
732  t_traction_cx, MultiPointRhsType::U);
733  for (int jj = 0; jj != 3; ++jj) {
734  auto v = t_rhs_tmp(jj).imag();
735  t_res_dP(jj, ii) = v / eps;
736  }
737  }
738  } else {
739 
741  t_cP(i, j) = (c * t_grad_sdf_v(i)) * t_grad_sdf_v(j);
742  t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
743  t_res_dP(i, j) = t_cQ(i, j);
744  }
745 
746  const double alpha = t_w / 2.;
747  size_t rr = 0;
748  for (; rr != nb_rows / 3; ++rr) {
749 
750  auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
751  auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
752 
753  for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
754  auto col_base = t_col_base(i) * t_material_normal(i);
755  const double beta = alpha * t_row_base * col_base;
756  t_mat(i, j) += beta * t_res_dP(i, j);
757  ++t_col_base;
758  ++t_mat;
759  }
760 
761  ++t_row_base;
762  }
763  for (; rr < nb_face_functions; ++rr)
764  ++t_row_base;
765 
766  next();
767  }
768  }
769 
771 }

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<ContactOps::CommonData> EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::commonDataPtr
private

Definition at line 172 of file EshelbianContact.hpp.

◆ contactTreePtr

boost::shared_ptr<ContactTree> EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::contactTreePtr
private

Definition at line 173 of file EshelbianContact.hpp.

◆ sdfMapRangePtr

boost::shared_ptr<std::map<int, Range> > EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::sdfMapRangePtr
private

Definition at line 174 of file EshelbianContact.hpp.


The documentation for this struct was generated from the following files:
EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::commonDataPtr
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
Definition: EshelbianContact.hpp:172
MoFEM::OpBaseImpl::locMat
MatrixDouble locMat
local entity block matrix
Definition: FormsIntegrators.hpp:239
FTensor::Tensor1< double, 3 >
ContactOps
Definition: EshelbianContact.hpp:10
ContactOps::constrain
double constrain(double sdf, double tn)
constrain function
Definition: ContactOps.hpp:572
FTensor::Kronecker_Delta
Kronecker Delta class.
Definition: Kronecker_Delta.hpp:15
EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::contactTreePtr
boost::shared_ptr< ContactTree > contactTreePtr
Definition: EshelbianContact.hpp:173
EshelbianPlasticity::checkSdf
auto checkSdf(EntityHandle fe_ent, std::map< int, Range > &sdf_map_range)
Definition: EshelbianContact.cpp:27
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
FTensor::Tensor2< double, 3, 3 >
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
ContactOps::AssemblyBoundaryEleOp
FormsIntegrators< BoundaryEleOp >::Assembly< A >::OpBase AssemblyBoundaryEleOp
Definition: EshelbianContact.hpp:17
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
FTensor::Index< 'i', 3 >
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
EshelbianPlasticity::OpConstrainBoundaryL2Lhs_dP::sdfMapRangePtr
boost::shared_ptr< std::map< int, Range > > sdfMapRangePtr
Definition: EshelbianContact.hpp:174
FTensor::kronecker_delta
Tensor2_Expr< Kronecker_Delta< T >, T, Dim0, Dim1, i, j > kronecker_delta(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
Rank 2.
Definition: Kronecker_Delta.hpp:81
EshelbianPlasticity::multiPointRhs
auto multiPointRhs(ContactTree::FaceData *face_data_ptr, FTensor::Tensor1< T1, 3 > &t_coords, FTensor::Tensor1< T2, 3 > &t_spatial_coords, FTensor::Tensor1< T3, 3 > &t_master_traction, MultiPointRhsType type, bool debug=false)
Definition: EshelbianContact.cpp:196
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
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
EshelbianPlasticity::U
@ U
Definition: EshelbianContact.cpp:193
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
EshelbianPlasticity::getSdf
auto getSdf(OP_PTR op_ptr, MatrixDouble &contact_disp, int block_id, bool eval_hessian)
Definition: EshelbianContact.cpp:36