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

Specialization for DisplacementCubitBcData. More...

#include <src/boundary_conditions/EssentialMPCsData.hpp>

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

Public Member Functions

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

Protected Attributes

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

Detailed Description

Specialization for DisplacementCubitBcData.

Template Parameters

Definition at line 81 of file EssentialMPCsData.hpp.

Constructor & Destructor Documentation

◆ EssentialPostProcRhs()

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

Definition at line 532 of file EssentialMPCsData.cpp.

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

◆ ~EssentialPostProcRhs()

Member Function Documentation

◆ operator()()

Definition at line 537 of file EssentialMPCsData.cpp.

537  {
538  MOFEM_LOG_CHANNEL("WORLD");
540 
541  if (auto fe_method_ptr = fePtr.lock()) {
542 
543  auto bc_mng = mField.getInterface<BcManager>();
544  auto is_mng = mField.getInterface<ISManager>();
545  auto vec_mng = mField.getInterface<VecManager>();
546  const auto problem_name = fe_method_ptr->problemPtr->getName();
547  // get connectivity from edge and assign master, slave ents
548  for (auto bc : bc_mng->getBcMapByBlockName()) {
549  if (auto mpc_bc = bc.second->mpcPtr) {
550 
551  auto &bc_id = bc.first;
552 
553  auto regex_str = (boost::format("%s_(.*)") % problem_name).str();
554  if (std::regex_match(bc_id, std::regex(regex_str))) {
555 
556  auto [field_name, block_name] =
557  BcManager::extractStringFromBlockId(bc_id, problem_name);
558 
559  auto get_field_coeffs = [&](auto field_name) {
560  auto field_ptr = mField.get_field_structure(field_name);
561  return field_ptr->getNbOfCoeffs();
562  };
563 
564  const auto nb_field_coeffs = get_field_coeffs(field_name);
565  constexpr auto max_nb_dofs_per_node = 6;
566 
567  if (nb_field_coeffs > max_nb_dofs_per_node)
568  MOFEM_LOG("WORLD", Sev::error)
569  << "MultiPointConstraints PreProcLhs<MPCsType>: support only "
570  "up to 6 dofs per node for now.";
571 
572  MOFEM_LOG("WORLD", Sev::noisy)
573  << "Apply MultiPointConstraints PreProc<MPCsType>: "
574  << problem_name << "_" << field_name << "_" << block_name;
575 
576  auto mpc_type = mpc_bc->mpcType;
577  switch (mpc_type) {
578  case MPC::COUPLING:
579  case MPC::TIE:
580  case MPC::RIGID_BODY:
581  break;
582  default:
583  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
584  "MPC type not implemented");
585  }
586 
587  auto master_verts = mpc_bc->mpcMasterEnts.subset_by_type(MBVERTEX);
588  auto slave_verts = mpc_bc->mpcSlaveEnts.subset_by_type(MBVERTEX);
589 
590  auto prb_name = fe_method_ptr->problemPtr->getName();
591 
592  auto get_flag = [&](int idx) { return (&mpc_bc->data.flag1)[idx]; };
593 
594  SmartPetscObj<IS> is_m[max_nb_dofs_per_node];
595  SmartPetscObj<IS> is_s[max_nb_dofs_per_node];
596 
597  for (int dd = 0; dd != nb_field_coeffs; dd++) {
598  if (get_flag(dd)) {
599 
600  CHKERR is_mng->isCreateProblemFieldAndRankLocal(
601  prb_name, ROW, field_name, dd, dd, is_m[dd], &master_verts);
602  CHKERR is_mng->isCreateProblemFieldAndRankLocal(
603  prb_name, ROW, field_name, dd, dd, is_s[dd], &slave_verts);
604 
605  }
606  }
607 
608  {
609 
610  if (auto fe_ptr = fePtr.lock()) {
611  auto snes_ctx = fe_ptr->snes_ctx;
612  auto ts_ctx = fe_ptr->ts_ctx;
613  const bool is_nonlinear = snes_ctx != FEMethod::CTX_SNESNONE ||
615  // is_nonlinear = is_nonlinear || mpc_bc->isReprocitical;
616  const int contrb = mpc_bc->isReprocitical ? 1 : 0;
617  SmartPetscObj<Vec> F =
618  vRhs ? vRhs : SmartPetscObj<Vec>(fe_ptr->f, true);
619 
620  if (fe_ptr->vecAssembleSwitch) {
621  if ((*fe_ptr->vecAssembleSwitch) && !vRhs) {
622  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
623  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
624  CHKERR VecAssemblyBegin(F);
625  CHKERR VecAssemblyEnd(F);
626  *fe_ptr->vecAssembleSwitch = false;
627  }
628  }
629 
630  auto vec_set_values = [&](auto is_xyz_m, auto is_xyz_s, double val) {
632  const int *m_index_ptr;
633  CHKERR ISGetIndices(is_xyz_m, &m_index_ptr);
634  const int *s_index_ptr;
635  CHKERR ISGetIndices(is_xyz_s, &s_index_ptr);
636  int size_m;
637  CHKERR ISGetLocalSize(is_xyz_m, &size_m);
638  int size_s;
639  CHKERR ISGetLocalSize(is_xyz_s, &size_s);
640 
641  double *f;
642  CHKERR VecGetArray(F, &f);
643  // if (size_m > size_s)
644  // MOFEM_LOG("WORLD", Sev::error)
645  // << "Size of master IS > Size of slave IS : " << size_m
646  // << " > " << size_s;
647  // if (size_m == 0)
648  // MOFEM_LOG("WORLD", Sev::error)
649  // << "Size of master IS is " << size_m;
650 
651  if (is_nonlinear) {
652  auto x = fe_ptr->x;
653  auto tmp_x = vectorDuplicate(F);
654 
655  CHKERR vec_mng->setLocalGhostVector(problem_name, ROW,
656  tmp_x, INSERT_VALUES,
657  SCATTER_FORWARD);
658  const double *v;
659  const double *u;
660 
661  CHKERR VecGetArrayRead(tmp_x, &u);
662  CHKERR VecGetArrayRead(x, &v);
663  if (size_m && size_s) {
664  for (auto i = 0; i != size_s; ++i) {
665  auto m_idx = size_m == 1 ? 0 : i;
666  f[s_index_ptr[i]] = val * (v[s_index_ptr[i]] - v[m_index_ptr[m_idx]]) + f[s_index_ptr[i]] * contrb ;
667  }
668  }
669  CHKERR VecRestoreArrayRead(x, &v);
670  CHKERR VecRestoreArrayRead(tmp_x, &u);
671  } else {
672  if (size_m && size_s)
673  for (auto i = 0; i != size_s; ++i) {
674  f[s_index_ptr[i]] = 0;
675  }
676  }
677 
678  CHKERR VecRestoreArray(F, &f);
679  CHKERR ISRestoreIndices(is_xyz_m, &m_index_ptr);
680  CHKERR ISRestoreIndices(is_xyz_s, &s_index_ptr);
681 
683  };
684 
685  for (int dd = 0; dd != nb_field_coeffs; dd++)
686  if (get_flag(dd)) {
687 
688  if (!mpc_bc->isReprocitical) {
689  CHKERR vec_set_values(is_m[dd], is_s[dd], vDiag);
690  } else {
691  auto &pn = mpc_bc->mPenalty;
692  CHKERR vec_set_values(is_m[dd], is_s[dd], pn);
693  CHKERR vec_set_values(is_s[dd], is_m[dd], pn);
694  }
695  }
696  };
697 
698  // User is responsible for assembly if vLhs is provided
699 
700  // ISView(is_xyz_m, PETSC_VIEWER_STDOUT_WORLD);
701  // CHKERR MatZeroRowsColumnsIS(B, is_xyz_m, vDiag, PETSC_NULL,
702  // PETSC_NULL);
703  }
704  }
705  }
706  }
707  } else {
708  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
709  "Can not lock shared pointer");
710  }
711 
713 }

Member Data Documentation

◆ fePtr

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

Definition at line 92 of file EssentialMPCsData.hpp.

◆ mField

Definition at line 91 of file EssentialMPCsData.hpp.

◆ vDiag

Definition at line 93 of file EssentialMPCsData.hpp.

◆ vRhs

Definition at line 94 of file EssentialMPCsData.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:460
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
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
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
MoFEM::EssentialPostProcRhs< MPCsType >::fePtr
boost::weak_ptr< FEMethod > fePtr
Definition: EssentialMPCsData.hpp:92
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
ROW
@ ROW
Definition: definitions.h:136
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::SnesMethod::CTX_SNESNONE
@ CTX_SNESNONE
Definition: LoopMethods.hpp:107
MoFEM::MPC::RIGID_BODY
@ RIGID_BODY
MoFEM::TSMethod::CTX_TSNONE
@ CTX_TSNONE
Definition: LoopMethods.hpp:145
MoFEM::EssentialPostProcRhs< MPCsType >::mField
MoFEM::Interface & mField
Definition: EssentialMPCsData.hpp:91
MoFEM::MPC::TIE
@ TIE
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
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
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:221
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::EssentialPostProcRhs< MPCsType >::vRhs
SmartPetscObj< Vec > vRhs
Definition: EssentialMPCsData.hpp:94
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::EssentialPostProcRhs< MPCsType >::vDiag
double vDiag
Definition: EssentialMPCsData.hpp:93
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::MPC::COUPLING
@ COUPLING
F
@ F
Definition: free_surface.cpp:394