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

Specialization for DisplacementCubitBcData. More...

#include <src/boundary_conditions/EssentialDisplacementCubitBcData.hpp>

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

Public Member Functions

 EssentialPostProcRhs (MoFEM::Interface &m_field, boost::shared_ptr< FEMethod > fe_ptr, double diag, SmartPetscObj< Vec > rhs=nullptr)
 
MoFEMErrorCode operator() ()
 

Protected Attributes

MoFEM::InterfacemField
 
boost::weak_ptr< FEMethodfePtr
 
double vDiag
 
SmartPetscObj< Vec > vRhs
 

Detailed Description

Specialization for DisplacementCubitBcData.

Template Parameters
Examples
adolc_plasticity.cpp, nonlinear_elastic.cpp, plastic.cpp, and seepage.cpp.

Definition at line 111 of file EssentialDisplacementCubitBcData.hpp.

Constructor & Destructor Documentation

◆ EssentialPostProcRhs()

MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >::EssentialPostProcRhs ( MoFEM::Interface m_field,
boost::shared_ptr< FEMethod fe_ptr,
double  diag,
SmartPetscObj< Vec >  rhs = nullptr 
)

Definition at line 171 of file EssentialDisplacementCubitBcData.cpp.

174  : mField(m_field), fePtr(fe_ptr), vDiag(diag), vRhs(rhs) {}

Member Function Documentation

◆ operator()()

Definition at line 176 of file EssentialDisplacementCubitBcData.cpp.

176  {
177  MOFEM_LOG_CHANNEL("WORLD");
179 
180  if (auto fe_method_ptr = fePtr.lock()) {
181 
182  auto bc_mng = mField.getInterface<BcManager>();
183  auto vec_mng = mField.getInterface<VecManager>();
184  auto is_mng = mField.getInterface<ISManager>();
185 
186  const auto problem_name = fe_method_ptr->problemPtr->getName();
187 
188  SmartPetscObj<IS> is_sum;
189 
190  for (auto bc : bc_mng->getBcMapByBlockName()) {
191  if (auto disp_bc = bc.second->dispBcPtr) {
192 
193  auto &bc_id = bc.first;
194 
195  auto regex_str = (boost::format("%s_(.*)") % problem_name).str();
196  if (std::regex_match(bc_id, std::regex(regex_str))) {
197 
198  auto [field_name, block_name] =
199  BcManager::extractStringFromBlockId(bc_id, problem_name);
200 
201  MOFEM_LOG("WORLD", Sev::noisy)
202  << "Apply EssentialPreProc<DisplacementCubitBcData>: "
203  << problem_name << "_" << field_name << "_" << block_name;
204 
205  const bool is_rotation =
206  disp_bc->data.flag4 || disp_bc->data.flag5 || disp_bc->data.flag6;
207 
208  auto ents = bc.second->bcEnts;
209 
210  std::array<SmartPetscObj<IS>, 3> is_xyz;
211  auto prb_name = fe_method_ptr->problemPtr->getName();
212 
213  if (disp_bc->data.flag1 || is_rotation) {
214  CHKERR is_mng->isCreateProblemFieldAndRankLocal(
215  prb_name, ROW, field_name, 0, 0, is_xyz[0], &ents);
216  }
217  if (disp_bc->data.flag2 || is_rotation) {
218  CHKERR is_mng->isCreateProblemFieldAndRankLocal(
219  prb_name, ROW, field_name, 1, 1, is_xyz[1], &ents);
220  }
221  if (disp_bc->data.flag3 || is_rotation) {
222  CHKERR is_mng->isCreateProblemFieldAndRankLocal(
223  prb_name, ROW, field_name, 2, 2, is_xyz[2], &ents);
224  }
225 
226  auto get_is_sum = [](auto is1, auto is2) {
227  IS is;
228  CHK_THROW_MESSAGE(ISExpand(is1, is2, &is), "is sum");
229  return SmartPetscObj<IS>(is);
230  };
231 
232  for (auto &is : is_xyz) {
233  if (is) {
234  if (!is_sum) {
235  is_sum = is;
236  } else {
237  is_sum = get_is_sum(is_sum, is);
238  }
239  }
240  }
241  }
242  }
243  }
244 
245  if (is_sum) {
246  if (auto fe_ptr = fePtr.lock()) {
247 
248  auto snes_ctx = fe_ptr->snes_ctx;
249  auto ts_ctx = fe_ptr->ts_ctx;
250 
251  SmartPetscObj<Vec> f =
252  vRhs ? vRhs : SmartPetscObj<Vec>(fe_ptr->f, true);
253 
254  if (fe_ptr->vecAssembleSwitch) {
255  if ((*fe_ptr->vecAssembleSwitch) && !vRhs) {
256  CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
257  CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
258  CHKERR VecAssemblyBegin(f);
259  CHKERR VecAssemblyEnd(f);
260  *fe_ptr->vecAssembleSwitch = false;
261  }
262  }
263 
264  const int *index_ptr;
265  CHKERR ISGetIndices(is_sum, &index_ptr);
266  int size;
267  CHKERR ISGetLocalSize(is_sum, &size);
268  double *a;
269  CHKERR VecGetArray(f, &a);
270 
271  auto tmp_x = vectorDuplicate(f);
272  CHKERR vec_mng->setLocalGhostVector(problem_name, ROW, tmp_x,
273  INSERT_VALUES, SCATTER_FORWARD);
274  const double *u;
275  CHKERR VecGetArrayRead(tmp_x, &u);
276 
277  if (snes_ctx != FEMethod::CTX_SNESNONE ||
279 
280  auto x = fe_ptr->x;
281  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
282  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
283 
284  const double *v;
285  CHKERR VecGetArrayRead(x, &v);
286 
287  for (auto i = 0; i != size; ++i) {
288  a[index_ptr[i]] = vDiag * (v[index_ptr[i]] - u[index_ptr[i]]);
289  }
290 
291  CHKERR VecRestoreArrayRead(x, &v);
292 
293  } else {
294  for (auto i = 0; i != size; ++i) {
295  a[index_ptr[i]] = vDiag * u[index_ptr[i]];
296  }
297  }
298 
299  CHKERR VecRestoreArrayRead(tmp_x, &u);
300  CHKERR VecRestoreArray(f, &a);
301  CHKERR ISRestoreIndices(is_sum, &index_ptr);
302  }
303  }
304 
305  } else {
306  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
307  "Can not lock shared pointer");
308  }
309 
311 }

Member Data Documentation

◆ fePtr

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

Definition at line 120 of file EssentialDisplacementCubitBcData.hpp.

◆ mField

Definition at line 119 of file EssentialDisplacementCubitBcData.hpp.

◆ vDiag

Definition at line 121 of file EssentialDisplacementCubitBcData.hpp.

◆ vRhs

Definition at line 122 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 refernce 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:1389
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:596
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
ROW
@ ROW
Definition: definitions.h:123
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::SnesMethod::CTX_SNESNONE
@ CTX_SNESNONE
Definition: LoopMethods.hpp:107
MoFEM::TSMethod::CTX_TSNONE
@ CTX_TSNONE
Definition: LoopMethods.hpp:145
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >::fePtr
boost::weak_ptr< FEMethod > fePtr
Definition: EssentialDisplacementCubitBcData.hpp:120
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >::vRhs
SmartPetscObj< Vec > vRhs
Definition: EssentialDisplacementCubitBcData.hpp:122
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
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >::vDiag
double vDiag
Definition: EssentialDisplacementCubitBcData.hpp:121
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
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:217
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >::mField
MoFEM::Interface & mField
Definition: EssentialDisplacementCubitBcData.hpp:119