18 #include "libmesh/quadrature.h" 19 #include "libmesh/elem.h" 20 #include "libmesh/point.h" 42 const Elem * primal_elem,
43 unsigned int sub_elem_index,
44 const QBase & qrule_msm,
45 std::vector<Point> & q_pts);
68 template <
typename Iterators,
typename Consumers,
typename ActionFunctor>
71 const Iterators & secondary_elems_to_mortar_segments,
77 const Consumers & consumers,
79 const std::map<
SubdomainID, std::deque<MaterialBase *>> & secondary_ip_sub_to_mats,
80 const std::map<
SubdomainID, std::deque<MaterialBase *>> & primary_ip_sub_to_mats,
81 const std::deque<MaterialBase *> & secondary_boundary_mats,
82 const ActionFunctor act,
83 const bool reinit_mortar_user_objects)
87 const auto primary_boundary_id = primary_secondary_boundary_id_pair.first;
88 const auto secondary_boundary_id = primary_secondary_boundary_id_pair.second;
91 unsigned int secondary_sub_elem_index = 0, primary_sub_elem_index = 0;
94 secondary_sub_elem_index = amg.
mortarSegmentMesh().get_elem_integer_index(
"secondary_sub_elem");
95 primary_sub_elem_index = amg.
mortarSegmentMesh().get_elem_integer_index(
"primary_sub_elem");
103 const auto & JxW_msm = assembly.
jxWMortar();
106 std::unordered_set<unsigned int> needed_mat_props;
107 for (
const auto & consumer : consumers)
109 const auto & mp_deps = consumer->getMatPropDependencies();
110 needed_mat_props.insert(mp_deps.begin(), mp_deps.end());
116 for (
const auto elem_to_msm : secondary_elems_to_mortar_segments)
118 const Elem * secondary_face_elem = subproblem.
mesh().
getMesh().elem_ptr(elem_to_msm->first);
120 const Elem * secondary_ip = secondary_face_elem->interior_parent();
121 unsigned int secondary_side_id = secondary_ip->which_side_am_i(secondary_face_elem);
122 const auto & secondary_ip_mats =
123 libmesh_map_find(secondary_ip_sub_to_mats, secondary_ip->subdomain_id());
125 const auto & msm_elems = elem_to_msm->second;
133 std::vector<Point> secondary_xi_pts, primary_xi_pts;
135 std::vector<Real> JxW;
138 unsigned int expected_length = 0;
142 for (
const auto msm_elem : msm_elems)
150 if (msm_elem->dim() == 1)
152 for (
unsigned int qp = 0; qp < qrule_msm->n_points(); qp++)
154 const Real eta = qrule_msm->qp(qp)(0);
157 const Real xi1_eta = 0.5 * (1 - eta) * msinfo.
xi1_a + 0.5 * (1 + eta) * msinfo.
xi1_b;
158 secondary_xi_pts.push_back(xi1_eta);
161 const Real xi2_eta = 0.5 * (1 - eta) * msinfo.
xi2_a + 0.5 * (1 + eta) * msinfo.
xi2_b;
162 primary_xi_pts.push_back(xi2_eta);
170 secondary_sub_elem_index,
174 msm_elem, msinfo.
primary_elem, primary_sub_elem_index, *qrule_msm, primary_xi_pts);
179 std::copy(std::begin(JxW_msm), std::end(JxW_msm), std::back_inserter(JxW));
183 expected_length += qrule_msm->n_points();
184 mooseAssert(secondary_xi_pts.size() == expected_length,
185 "Fewer than expected secondary quadrature points");
186 mooseAssert(primary_xi_pts.size() == expected_length,
187 "Fewer than expected primary quadrature points");
190 mooseAssert(JxW.size() == expected_length,
"Fewer than expected JxW values computed");
197 assembly.
reinitDual(secondary_face_elem, secondary_xi_pts, JxW);
199 unsigned int n_segment = 0;
202 for (
const auto msm_elem : msm_elems)
207 std::vector<Point> xi1_pts, xi2_pts;
213 const Elem * primary_ip = msinfo.
primary_elem->interior_parent();
214 unsigned int primary_side_id = primary_ip->which_side_am_i(msinfo.
primary_elem);
215 const auto & primary_ip_mats =
216 libmesh_map_find(primary_ip_sub_to_mats, primary_ip->subdomain_id());
224 const unsigned int start = (n_segment - 1) * qrule_msm->n_points();
225 const unsigned int end = n_segment * qrule_msm->n_points();
227 xi1_pts.begin(), secondary_xi_pts.begin() + start, secondary_xi_pts.begin() + end);
228 xi2_pts.insert(xi2_pts.begin(), primary_xi_pts.begin() + start, primary_xi_pts.begin() + end);
230 const Elem * reinit_secondary_elem = secondary_ip;
235 reinit_secondary_elem = fe_problem.
mesh().
elemPtr(reinit_secondary_elem->id());
251 reinit_secondary_elem, secondary_side_id, TOLERANCE, &xi1_pts,
nullptr, tid);
253 const Elem * reinit_primary_elem = primary_ip;
258 reinit_primary_elem = fe_problem.
mesh().
elemPtr(reinit_primary_elem->id());
262 reinit_primary_elem, primary_side_id, TOLERANCE, &xi2_pts,
nullptr, tid);
290 secondary_boundary_id, tid,
false, &secondary_boundary_mats);
292 if (reinit_mortar_user_objects)
315 template <
typename Consumers>
321 std::map<
SubdomainID, std::deque<MaterialBase *>> & secondary_ip_sub_to_mats,
322 std::map<
SubdomainID, std::deque<MaterialBase *>> & primary_ip_sub_to_mats,
323 std::deque<MaterialBase *> & secondary_boundary_mats)
325 secondary_ip_sub_to_mats.clear();
326 primary_ip_sub_to_mats.clear();
327 secondary_boundary_mats.clear();
330 auto get_required_sub_mats =
331 [&mat_warehouse, tid, &consumers](
335 if (mat_warehouse[mat_data_type].hasActiveBlockObjects(sub_id, tid))
346 for (
const auto secondary_ip_sub : secondary_ip_sub_ids)
347 secondary_ip_sub_to_mats.emplace(
352 for (
const auto primary_ip_sub : primary_ip_sub_ids)
353 primary_ip_sub_to_mats.emplace(
358 const auto secondary_boundary = boundary_pr.second;
359 if (mat_warehouse.hasActiveBoundaryObjects(secondary_boundary, tid))
361 auto & boundary_mats = mat_warehouse.getActiveBoundaryObjects(secondary_boundary, tid);
362 secondary_boundary_mats =
virtual MooseMesh & mesh()=0
void setActiveMaterialProperties(const std::unordered_set< unsigned int > &mat_prop_ids, const THREAD_ID tid)
Record and set the material properties required by the current computing thread.
void reinitNeighborLowerDElem(const Elem *elem, const THREAD_ID tid=0)
reinitialize a neighboring lower dimensional element
const std::pair< BoundaryID, BoundaryID > & primarySecondaryBoundaryIDPair() const
virtual void reinitLowerDElem(const Elem *lower_d_elem, const THREAD_ID tid, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)
Keeps track of stuff related to assembling.
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
virtual Elem * elemPtr(const dof_id_type i)
void setupMortarMaterials(const Consumers &consumers, FEProblemBase &fe_problem, const AutomaticMortarGeneration &amg, const THREAD_ID tid, std::map< SubdomainID, std::deque< MaterialBase *>> &secondary_ip_sub_to_mats, std::map< SubdomainID, std::deque< MaterialBase *>> &primary_ip_sub_to_mats, std::deque< MaterialBase *> &secondary_boundary_mats)
This function creates containers of materials necessary to execute the mortar method for a supplied s...
MaterialDataType
MaterialData types.
void reinitMortarUserObjects(BoundaryID primary_boundary_id, BoundaryID secondary_boundary_id, bool displaced)
Call reinit on mortar user objects with matching primary boundary ID, secondary boundary ID...
void reinitMaterialsBoundary(BoundaryID boundary_id, const THREAD_ID tid, bool swap_stateful=true, const std::deque< MaterialBase *> *reinit_mats=nullptr)
reinit materials on a boundary
const libMesh::QBase *const & qRuleMortar() const
Returns a reference to the quadrature rule for the mortar segments.
const std::set< SubdomainID > & primaryIPSubIDs() const
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
const MaterialWarehouse & getRegularMaterialsWarehouse() const
void reinitMaterialsFace(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true, const std::deque< MaterialBase *> *reinit_mats=nullptr)
reinit materials on element faces
const std::set< SubdomainID > & secondaryIPSubIDs() const
This class is a container/interface for the objects involved in automatic generation of mortar spaces...
virtual void reinitElemFaceRef(const Elem *elem, unsigned int side, Real tolerance, const std::vector< Point > *const pts, const std::vector< Real > *const weights=nullptr, const THREAD_ID tid=0) override
reinitialize FE objects on a given element on a given side at a given set of reference points and the...
const Elem * primary_elem
bool needDual() const
Indicates whether dual shape functions are used (computation is now repeated on each element so expen...
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
void reinitMaterialsNeighbor(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true, const std::deque< MaterialBase *> *reinit_mats=nullptr)
reinit materials on the neighboring element face
virtual void reinitNeighborFaceRef(const Elem *neighbor_elem, unsigned int neighbor_side, Real tolerance, const std::vector< Point > *const pts, const std::vector< Real > *const weights=nullptr, const THREAD_ID tid=0) override
reinitialize FE objects on a given neighbor element on a given side at a given set of reference point...
void loopOverMortarSegments(const Iterators &secondary_elems_to_mortar_segments, Assembly &assembly, SubProblem &subproblem, FEProblemBase &fe_problem, const AutomaticMortarGeneration &amg, const bool displaced, const Consumers &consumers, const THREAD_ID tid, const std::map< SubdomainID, std::deque< MaterialBase *>> &secondary_ip_sub_to_mats, const std::map< SubdomainID, std::deque< MaterialBase *>> &primary_ip_sub_to_mats, const std::deque< MaterialBase *> &secondary_boundary_mats, const ActionFunctor act, const bool reinit_mortar_user_objects)
This method will loop over pairs of secondary elements and their corresponding mortar segments...
const Elem * secondary_elem
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Holds xi^(1), xi^(2), and other data for a given mortar segment.
Generic class for solving transient nonlinear problems.
static std::deque< MaterialBase * > buildRequiredMaterials(const Consumers &mat_consumers, const std::vector< std::shared_ptr< MaterialBase >> &mats, const bool allow_stateful)
Build the materials required by a set of consumer objects.
const std::unordered_map< const Elem *, MortarSegmentInfo > & mortarSegmentMeshElemToInfo() const
virtual MooseMesh & mesh() override
const MeshBase & mortarSegmentMesh() const
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
const std::vector< Real > & jxWMortar() const
Returns a reference to JxW for mortar segment elements.
void reinitMortarElem(const Elem *elem, const THREAD_ID tid=0)
Reinit a mortar element to obtain a valid JxW.
void projectQPoints3d(const Elem *msm_elem, const Elem *primal_elem, unsigned int sub_elem_index, const QBase &qrule_msm, std::vector< Point > &q_pts)
3D projection operator for mapping qpoints on mortar segments to secondary or primary elements ...
void reinitDual(const Elem *elem, const std::vector< Point > &pts, const std::vector< Real > &JxW)
Reintialize dual basis coefficients based on a customized quadrature rule.