v0.14.0
Public Member Functions | Protected Attributes | List of all members
MoFEM::EssentialPreProc< DisplacementCubitBcData > Struct Reference

Specialization for DisplacementCubitBcData. More...

#include <src/boundary_conditions/EssentialDisplacementCubitBcData.hpp>

Collaboration diagram for MoFEM::EssentialPreProc< DisplacementCubitBcData >:
[legend]

Public Member Functions

 EssentialPreProc (MoFEM::Interface &m_field, boost::shared_ptr< FEMethod > fe_ptr, std::vector< boost::shared_ptr< ScalingMethod >> smv, bool get_coords=false)
 
virtual ~EssentialPreProc ()=default
 
MoFEMErrorCode operator() ()
 

Protected Attributes

MoFEM::InterfacemField
 
boost::weak_ptr< FEMethodfePtr
 
VecOfTimeScalingMethods vecOfTimeScalingMethods
 
bool getCoords
 

Detailed Description

Specialization for DisplacementCubitBcData.

Specialization to enforce blocksets which DisplacementCubitBcData ptr. That is to enforce displacement constraints. set

Template Parameters
Examples
adolc_plasticity.cpp, dynamic_first_order_con_law.cpp, plastic.cpp, seepage.cpp, and thermo_elastic.cpp.

Definition at line 91 of file EssentialDisplacementCubitBcData.hpp.

Constructor & Destructor Documentation

◆ EssentialPreProc()

MoFEM::EssentialPreProc< DisplacementCubitBcData >::EssentialPreProc ( MoFEM::Interface m_field,
boost::shared_ptr< FEMethod fe_ptr,
std::vector< boost::shared_ptr< ScalingMethod >>  smv,
bool  get_coords = false 
)

Definition at line 15 of file EssentialDisplacementCubitBcData.cpp.

18  : mField(m_field), fePtr(fe_ptr), vecOfTimeScalingMethods(smv),
19  getCoords(get_coords) {}

◆ ~EssentialPreProc()

Member Function Documentation

◆ operator()()

Definition at line 21 of file EssentialDisplacementCubitBcData.cpp.

21  {
22  MOFEM_LOG_CHANNEL("WORLD");
24 
25  if (auto fe_method_ptr = fePtr.lock()) {
26 
27  auto bc_mng = mField.getInterface<BcManager>();
28  auto fb = mField.getInterface<FieldBlas>();
29  const auto problem_name = fe_method_ptr->problemPtr->getName();
30 
31  for (auto bc : bc_mng->getBcMapByBlockName()) {
32  if (auto disp_bc = bc.second->dispBcPtr) {
33 
34  auto &bc_id = bc.first;
35 
36  auto regex_str = (boost::format("%s_(.*)") % problem_name).str();
37  if (std::regex_match(bc_id, std::regex(regex_str))) {
38 
39  auto [field_name, block_name] =
40  BcManager::extractStringFromBlockId(bc_id, problem_name);
41 
42  auto get_field_coeffs = [&](auto field_name) {
43  auto field_ptr = mField.get_field_structure(field_name);
44  return field_ptr->getNbOfCoeffs();
45  };
46  const auto nb_field_coeffs = get_field_coeffs(field_name);
47 
48  MOFEM_LOG("WORLD", Sev::noisy)
49  << "Apply EssentialPreProc<DisplacementCubitBcData>: "
50  << problem_name << "_" << field_name << "_" << block_name;
51 
52  FTensor::Tensor1<double, 3> t_angles{0., 0., 0.};
53  FTensor::Tensor1<double, 3> t_vals{0., 0., 0.};
54  FTensor::Tensor1<double, 3> t_off{0., 0., 0.};
55 
56  if (auto ext_disp_bc =
57  dynamic_cast<DisplacementCubitBcDataWithRotation const *>(
58  disp_bc.get())) {
59  for (int a = 0; a != 3; ++a)
60  t_off(a) = ext_disp_bc->rotOffset[a];
61  }
62 
63  auto scale_value = [&](const double &c) {
64  double val = c;
65  for (auto s : vecOfTimeScalingMethods) {
66  val *= s->getScale(fe_method_ptr->ts_t);
67  }
68  return val;
69  };
70 
71  if (disp_bc->data.flag1 == 1)
72  t_vals(0) = scale_value(-disp_bc->data.value1);
73  if (disp_bc->data.flag2 == 1)
74  t_vals(1) = scale_value(-disp_bc->data.value2);
75  if (disp_bc->data.flag3 == 1)
76  t_vals(2) = scale_value(-disp_bc->data.value3);
77  if (disp_bc->data.flag4 == 1)
78  t_angles(0) = scale_value(-disp_bc->data.value4);
79  if (disp_bc->data.flag5 == 1)
80  t_angles(1) = scale_value(-disp_bc->data.value5);
81  if (disp_bc->data.flag6 == 1)
82  t_angles(2) = scale_value(-disp_bc->data.value6);
83 
84  int coeff;
85  std::array<std::vector<double>, 3> coords;
86  int idx;
87 
88  const bool is_rotation =
89  disp_bc->data.flag4 || disp_bc->data.flag5 || disp_bc->data.flag6;
90 
91  auto lambda = [&](boost::shared_ptr<FieldEntity> field_entity_ptr) {
93 
94  auto v = t_vals(coeff);
95  if (is_rotation) {
97  coords[0][idx], coords[1][idx], coords[2][idx]);
99  t_angles, t_coords, t_off)(coeff);
100  }
101  if (getCoords) {
102  v += coords[coeff][idx];
103  }
104 
105  field_entity_ptr->getEntFieldData()[coeff] = v;
106  ++idx;
107 
109  };
110 
111  auto zero_lambda =
112  [&](boost::shared_ptr<FieldEntity> field_entity_ptr) {
114  auto size = field_entity_ptr->getEntFieldData().size();
115  for (int i = coeff; i < size; i += nb_field_coeffs)
116  field_entity_ptr->getEntFieldData()[i] = 0;
118  };
119 
120  auto verts = bc.second->bcEnts.subset_by_type(MBVERTEX);
121  auto not_verts = subtract(bc.second->bcEnts, verts);
122 
123  if (getCoords || is_rotation) {
124  for (auto d : {0, 1, 2})
125  coords[d].resize(verts.size());
126  CHKERR mField.get_moab().get_coords(verts, &*coords[0].begin(),
127  &*coords[1].begin(),
128  &*coords[2].begin());
129  }
130 
131  if (disp_bc->data.flag1 || disp_bc->data.flag5 ||
132  disp_bc->data.flag6) {
133  idx = 0;
134  coeff = 0;
135  CHKERR fb->fieldLambdaOnEntities(lambda, field_name, &verts);
136  CHKERR fb->fieldLambdaOnEntities(zero_lambda, field_name,
137  &not_verts);
138  }
139  if (disp_bc->data.flag2 || disp_bc->data.flag4 ||
140  disp_bc->data.flag6) {
141  idx = 0;
142  coeff = 1;
143  if (nb_field_coeffs > 1) {
144  CHKERR fb->fieldLambdaOnEntities(lambda, field_name, &verts);
145  CHKERR fb->fieldLambdaOnEntities(zero_lambda, field_name,
146  &not_verts);
147  }
148  }
149  if (disp_bc->data.flag3 || disp_bc->data.flag4 ||
150  disp_bc->data.flag5 || is_rotation) {
151  idx = 0;
152  coeff = 2;
153  if (nb_field_coeffs > 2) {
154  CHKERR fb->fieldLambdaOnEntities(lambda, field_name, &verts);
155  CHKERR fb->fieldLambdaOnEntities(zero_lambda, field_name,
156  &not_verts);
157  }
158  }
159  }
160  }
161  }
162 
163  } else {
164  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
165  "Can not lock shared pointer");
166  }
167 
169 }

Member Data Documentation

◆ fePtr

boost::weak_ptr<FEMethod> MoFEM::EssentialPreProc< DisplacementCubitBcData >::fePtr
protected

Definition at line 103 of file EssentialDisplacementCubitBcData.hpp.

◆ getCoords

Definition at line 105 of file EssentialDisplacementCubitBcData.hpp.

◆ mField

Definition at line 102 of file EssentialDisplacementCubitBcData.hpp.

◆ vecOfTimeScalingMethods

Definition at line 104 of file EssentialDisplacementCubitBcData.hpp.


The documentation for this struct was generated from the following files:
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::BcManager::extractStringFromBlockId
static std::pair< std::string, std::string > extractStringFromBlockId(const std::string block_id, const std::string prb_name)
Extract block name and block name form block id.
Definition: BcManager.cpp:1381
MoFEM::Field::getNbOfCoeffs
FieldCoefficientsNumber getNbOfCoeffs() const
Get number of field coefficients.
Definition: FieldMultiIndices.hpp:202
FTensor::Tensor1< double, 3 >
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
MoFEM::EssentialPreProc< DisplacementCubitBcData >::vecOfTimeScalingMethods
VecOfTimeScalingMethods vecOfTimeScalingMethods
Definition: EssentialDisplacementCubitBcData.hpp:104
MoFEM::CoreInterface::get_field_structure
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::EssentialPreProc< DisplacementCubitBcData >::mField
MoFEM::Interface & mField
Definition: EssentialDisplacementCubitBcData.hpp:102
MoFEM::EssentialPreProc< DisplacementCubitBcData >::getCoords
bool getCoords
Definition: EssentialDisplacementCubitBcData.hpp:105
MoFEM::DisplacementCubitBcDataWithRotation::GetRotDisp
static FTensor::Tensor1< double, 3 > GetRotDisp(const FTensor::Tensor1< double, 3 > &angles, FTensor::Tensor1< double, 3 > coordinates, FTensor::Tensor1< double, 3 > offset=FTensor::Tensor1< double, 3 >{ 0., 0., 0.})
Calculates the rotated displacement given the rotation angles, coordinates, and an optional offset.
Definition: EssentialDisplacementCubitBcData.hpp:40
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
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::EssentialPreProc< DisplacementCubitBcData >::fePtr
boost::weak_ptr< FEMethod > fePtr
Definition: EssentialDisplacementCubitBcData.hpp:103