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

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

Inheritance diagram for FractureMechanics::OpAnalyticalSpatialTraction:
[legend]
Collaboration diagram for FractureMechanics::OpAnalyticalSpatialTraction:
[legend]

Public Member Functions

 OpAnalyticalSpatialTraction (MoFEM::Interface &m_field, const bool &set_singular_coordinates, const Range &crack_front_nodes, const Range &crack_front_nodes_edges, const Range &crack_front_elements, PetscBool add_singularity, const std::string field_name, const Range tris, boost::ptr_vector< MethodForForceScaling > &methods_op, boost::ptr_vector< NeumannForcesSurface::MethodForAnalyticalForce > &analytical_force_op)
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 

Public Attributes

const Range tRis
 
boost::ptr_vector< MethodForForceScaling > & methodsOp
 
boost::ptr_vector< NeumannForcesSurface::MethodForAnalyticalForce > & analyticalForceOp
 
CrackFrontSingularBase< VolumeElementForcesAndSourcesCoreOnSide, VolumeElementForcesAndSourcesCorevolSideFe
 
VectorDouble sNrm
 Length of the normal vector. More...
 
bool setSingularCoordinates
 
VectorDouble Nf
 

Detailed Description

Definition at line 404 of file CrackFrontElement.hpp.

Constructor & Destructor Documentation

◆ OpAnalyticalSpatialTraction()

FractureMechanics::OpAnalyticalSpatialTraction::OpAnalyticalSpatialTraction ( MoFEM::Interface m_field,
const bool set_singular_coordinates,
const Range crack_front_nodes,
const Range crack_front_nodes_edges,
const Range crack_front_elements,
PetscBool  add_singularity,
const std::string  field_name,
const Range  tris,
boost::ptr_vector< MethodForForceScaling > &  methods_op,
boost::ptr_vector< NeumannForcesSurface::MethodForAnalyticalForce > &  analytical_force_op 
)

Definition at line 881 of file CrackFrontElement.cpp.

891  tRis(tris), methodsOp(methods_op), analyticalForceOp(analytical_force_op),
892  volSideFe(m_field, set_singular_coordinates, crack_front_nodes,
893  crack_front_nodes_edges, crack_front_elements, add_singularity),
894  setSingularCoordinates(set_singular_coordinates) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode FractureMechanics::OpAnalyticalSpatialTraction::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)

Definition at line 897 of file CrackFrontElement.cpp.

898  {
900 
901  if (data.getIndices().size() == 0)
903  EntityHandle ent = getNumeredEntFiniteElementPtr()->getEnt();
904  if (tRis.find(ent) == tRis.end())
906 
907  if (type == MBVERTEX) {
908  CHKERR loopSideVolumes("ELASTIC", volSideFe);
909  }
910 
911  int rank = data.getFieldDofs()[0]->getNbOfCoeffs();
912  int nb_row_dofs = data.getIndices().size() / rank;
913 
914  Nf.resize(data.getIndices().size(), false);
915  Nf.clear();
916 
917  VectorDouble3 coords(3);
918  VectorDouble3 normal(3);
919  VectorDouble3 force(3);
920  VectorDouble3 FTforce(3);
921 
925 
926  // If singular element apply transformation to h tensor
927  if (type == MBVERTEX && volSideFe.singularElement) {
928  double *inv_jac_ptr = &*volSideFe.invSJac.data().begin();
930  inv_jac_ptr, &inv_jac_ptr[1], &inv_jac_ptr[2], &inv_jac_ptr[3],
931  &inv_jac_ptr[4], &inv_jac_ptr[5], &inv_jac_ptr[6], &inv_jac_ptr[7],
932  &inv_jac_ptr[8], 9);
933  FTensor::Tensor1<double *, 3> t_normal = getFTensor1Normal();
934  double nrm = 2 * getArea();
935  FTensor::Tensor1<double, 3> t_s_normal;
936  sNrm.resize(data.getN().size1(), false);
937  for (unsigned int gg = 0; gg != data.getN().size1(); gg++) {
938  // Push normal to tranformed singular element
939  double det = volSideFe.detS[gg];
940  t_s_normal(i) = det * t_inv_s_jac(j, i) * t_normal(j) / nrm;
941  sNrm[gg] = sqrt(t_s_normal(i) * t_s_normal(i));
942  ++t_inv_s_jac;
943  }
944  }
945 
946  FTensor::Tensor1<double *, 3> t_force(&force[0], &force[1], &force[2]);
947 
948  const auto &coords_at_gauss_pts = getCoordsAtGaussPts();
949  const auto &normals_at_gauss_pts = getNormalsAtGaussPts();
950 
951  for (unsigned int gg = 0; gg < data.getN().size1(); gg++) {
952 
953  // get integration weight and Jacobian of integration point (area of face)
954  double val = getGaussPts()(2, gg);
955  val *= 0.5 * cblas_dnrm2(3, &getNormalsAtGaussPts()(gg, 0), 1);
956  for (int dd = 0; dd != 3; dd++) {
957  coords[dd] = coords_at_gauss_pts(gg, dd);
958  normal[dd] = normals_at_gauss_pts(gg, dd);
959  }
961  for (int dd = 0; dd != 3; dd++) {
962  coords[dd] += volSideFe.singularDisp(gg, dd);
963  }
964  val *= sNrm[gg];
965  }
966 
967  for (boost::ptr_vector<
969  analyticalForceOp.begin();
970  vit != analyticalForceOp.end(); vit++) {
971 
972  // Calculate force
973  CHKERR vit->getForce(ent, coords, normal, force);
974 
975  for (int rr = 0; rr != 3; rr++) {
976  cblas_daxpy(nb_row_dofs, val * force[rr], &data.getN()(gg, 0), 1,
977  &Nf[rr], rank);
978  }
979  }
980  }
981 
982  // Scale force using user defined scaling operator
984 
985  // Assemble force into vector
986  CHKERR VecSetValues(getFEMethod()->snes_f, data.getIndices().size(),
987  &data.getIndices()[0], &Nf[0], ADD_VALUES);
988 
990 }

Member Data Documentation

◆ analyticalForceOp

boost::ptr_vector<NeumannForcesSurface::MethodForAnalyticalForce>& FractureMechanics::OpAnalyticalSpatialTraction::analyticalForceOp

Definition at line 410 of file CrackFrontElement.hpp.

◆ methodsOp

boost::ptr_vector<MethodForForceScaling>& FractureMechanics::OpAnalyticalSpatialTraction::methodsOp

Definition at line 408 of file CrackFrontElement.hpp.

◆ Nf

VectorDouble FractureMechanics::OpAnalyticalSpatialTraction::Nf

Definition at line 429 of file CrackFrontElement.hpp.

◆ setSingularCoordinates

bool FractureMechanics::OpAnalyticalSpatialTraction::setSingularCoordinates

Definition at line 418 of file CrackFrontElement.hpp.

◆ sNrm

VectorDouble FractureMechanics::OpAnalyticalSpatialTraction::sNrm

Length of the normal vector.

Definition at line 417 of file CrackFrontElement.hpp.

◆ tRis

const Range FractureMechanics::OpAnalyticalSpatialTraction::tRis

Definition at line 407 of file CrackFrontElement.hpp.

◆ volSideFe

CrackFrontSingularBase<VolumeElementForcesAndSourcesCoreOnSide, VolumeElementForcesAndSourcesCore> FractureMechanics::OpAnalyticalSpatialTraction::volSideFe

Definition at line 414 of file CrackFrontElement.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
FractureMechanics::OpAnalyticalSpatialTraction::sNrm
VectorDouble sNrm
Length of the normal vector.
Definition: CrackFrontElement.hpp:417
MoFEM::Types::VectorDouble3
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
FractureMechanics::CrackFrontSingularBase::singularElement
bool singularElement
Definition: CrackFrontElement.hpp:47
FractureMechanics::OpAnalyticalSpatialTraction::analyticalForceOp
boost::ptr_vector< NeumannForcesSurface::MethodForAnalyticalForce > & analyticalForceOp
Definition: CrackFrontElement.hpp:410
FractureMechanics::OpAnalyticalSpatialTraction::methodsOp
boost::ptr_vector< MethodForForceScaling > & methodsOp
Definition: CrackFrontElement.hpp:408
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1576
FTensor::Tensor2< double *, 3, 3 >
FractureMechanics::CrackFrontSingularBase::invSJac
MatrixDouble invSJac
Definition: CrackFrontElement.hpp:42
NeumannForcesSurface::MethodForAnalyticalForce
Analytical force method.
Definition: SurfacePressure.hpp:21
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
convert.type
type
Definition: convert.py:64
FractureMechanics::CrackFrontSingularBase::singularDisp
MatrixDouble singularDisp
Definition: CrackFrontElement.hpp:43
FractureMechanics::OpAnalyticalSpatialTraction::volSideFe
CrackFrontSingularBase< VolumeElementForcesAndSourcesCoreOnSide, VolumeElementForcesAndSourcesCore > volSideFe
Definition: CrackFrontElement.hpp:414
FractureMechanics::CrackFrontSingularBase::detS
VectorDouble detS
Definition: CrackFrontElement.hpp:45
FractureMechanics::OpAnalyticalSpatialTraction::setSingularCoordinates
bool setSingularCoordinates
Definition: CrackFrontElement.hpp:418
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
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::OpAnalyticalSpatialTraction::tRis
const Range tRis
Definition: CrackFrontElement.hpp:407
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
FractureMechanics::OpAnalyticalSpatialTraction::Nf
VectorDouble Nf
Definition: CrackFrontElement.hpp:429
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MethodForForceScaling::applyScale
static MoFEMErrorCode applyScale(const FEMethod *fe, boost::ptr_vector< MethodForForceScaling > &methods_op, VectorDouble &nf)
Definition: MethodForForceScaling.hpp:21
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567