https://mooseframework.inl.gov
AutomaticMortarGeneration.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "MortarSegmentInfo.h"
13 #include "MooseHashing.h"
14 #include "ConsoleStreamInterface.h"
15 #include "MooseError.h"
16 #include "MooseUtils.h"
17 
18 // libMesh includes
19 #include "libmesh/id_types.h"
20 #include "libmesh/equation_systems.h"
21 #include "libmesh/elem.h"
22 #include "libmesh/int_range.h"
23 
24 // C++ includes
25 #include <set>
26 #include <memory>
27 #include <vector>
28 #include <unordered_map>
29 
30 // Forward declarations
31 namespace libMesh
32 {
33 class MeshBase;
34 class System;
35 }
36 class GetPot;
37 
38 // Using statements
42 using libMesh::Elem;
43 using libMesh::MeshBase;
44 using libMesh::Node;
45 using libMesh::Point;
46 using libMesh::Real;
48 
51 
57 {
58 public:
63  const static std::string system_name;
64 
70  MeshBase & mesh_in,
71  const std::pair<BoundaryID, BoundaryID> & boundary_key,
72  const std::pair<SubdomainID, SubdomainID> & subdomain_key,
73  bool on_displaced,
74  bool periodic,
75  const bool debug,
76  const bool correct_edge_dropping,
77  const Real minimum_projection_angle);
78 
85  void buildNodeToElemMaps();
86 
94  void computeNodalGeometry();
95 
109  void projectSecondaryNodes();
110 
125  void projectPrimaryNodes();
126 
143  void buildMortarSegmentMesh();
144 
157 
161  void msmStatistics();
162 
166  void clear();
167 
171  bool onDisplaced() const { return _on_displaced; }
172 
176  std::vector<Point> getNodalNormals(const Elem & secondary_elem) const;
177 
182  std::array<MooseUtils::SemidynamicVector<Point, 9>, 2>
183  getNodalTangents(const Elem & secondary_elem) const;
184 
189  std::map<unsigned int, unsigned int>
190  getSecondaryIpToLowerElementMap(const Elem & lower_secondary_elem) const;
191 
198  std::map<unsigned int, unsigned int>
199  getPrimaryIpToLowerElementMap(const Elem & primary_elem,
200  const Elem & primary_elem_ip,
201  const Elem & lower_secondary_elem) const;
202 
211  std::vector<Point> getNormals(const Elem & secondary_elem,
212  const std::vector<Point> & xi1_pts) const;
213 
222  std::vector<Point> getNormals(const Elem & secondary_elem,
223  const std::vector<Real> & oned_xi1_pts) const;
224 
232  const Elem * getSecondaryLowerdElemFromSecondaryElem(dof_id_type secondary_elem_id) const;
233 
238  void computeInactiveLMNodes();
239 
245 
250  void computeInactiveLMElems();
251 
255  const std::unordered_map<dof_id_type, std::unordered_set<dof_id_type>> &
257  {
259  }
260 
264  const std::pair<BoundaryID, BoundaryID> & primarySecondaryBoundaryIDPair() const;
265 
269  const MeshBase & mortarSegmentMesh() const { return *_mortar_segment_mesh; }
270 
274  const std::unordered_map<const Elem *, MortarSegmentInfo> & mortarSegmentMeshElemToInfo() const
275  {
276  return _msm_elem_to_info;
277  }
278 
279  int dim() const { return _mesh.mesh_dimension(); }
280 
284  const std::unordered_set<const Node *> & getInactiveLMNodes() const
285  {
287  }
288 
292  const std::unordered_set<const Elem *> & getInactiveLMElems() const
293  {
295  }
296 
298 
299  using MortarFilterIter =
300  std::unordered_map<dof_id_type, std::set<Elem *, CompareDofObjectsByID>>::const_iterator;
301 
308  std::vector<MortarFilterIter> secondariesToMortarSegments(const Node & node) const;
309 
314  const std::unordered_map<dof_id_type, std::set<Elem *, CompareDofObjectsByID>> &
316  {
318  }
319 
323  const std::set<SubdomainID> & secondaryIPSubIDs() const { return _secondary_ip_sub_ids; }
324 
328  const std::set<SubdomainID> & primaryIPSubIDs() const { return _primary_ip_sub_ids; }
329 
333  const std::unordered_map<dof_id_type, std::vector<const Elem *>> & nodesToSecondaryElem() const
334  {
336  }
337 
341  void initOutput();
342 
343 private:
348 
351 
354 
356  std::set<BoundaryID> _secondary_requested_boundary_ids;
357 
359  std::set<BoundaryID> _primary_requested_boundary_ids;
360 
363  std::vector<std::pair<BoundaryID, BoundaryID>> _primary_secondary_boundary_id_pairs;
364 
366  std::unordered_map<dof_id_type, std::vector<const Elem *>> _nodes_to_secondary_elem_map;
367  std::unordered_map<dof_id_type, std::vector<const Elem *>> _nodes_to_primary_elem_map;
368 
386  std::unordered_map<std::pair<const Node *, const Elem *>, std::pair<Real, const Elem *>>
388 
394  std::map<std::tuple<dof_id_type, const Node *, const Elem *>, std::pair<Real, const Elem *>>
396 
399  std::unique_ptr<MeshBase> _mortar_segment_mesh;
400 
404  std::unordered_map<const Elem *, MortarSegmentInfo> _msm_elem_to_info;
405 
408  std::unordered_map<const Elem *, unsigned int> _lower_elem_to_side_id;
409 
412  std::vector<std::pair<SubdomainID, SubdomainID>> _primary_secondary_subdomain_id_pairs;
413 
416  std::set<SubdomainID> _secondary_boundary_subdomain_ids;
417  std::set<SubdomainID> _primary_boundary_subdomain_ids;
418 
427  std::unordered_map<dof_id_type, std::unordered_set<dof_id_type>> _mortar_interface_coupling;
428 
430  std::unordered_map<const Node *, Point> _secondary_node_to_nodal_normal;
431 
434  std::unordered_map<const Node *, std::array<Point, 2>> _secondary_node_to_hh_nodal_tangents;
435 
437  std::unordered_map<dof_id_type, const Elem *> _secondary_element_to_secondary_lowerd_element;
438 
439  // List of inactive lagrange multiplier nodes (for nodal variables)
440  std::unordered_set<const Node *> _inactive_local_lm_nodes;
441 
443  std::unordered_set<const Elem *> _inactive_local_lm_elems;
444 
448  std::unordered_map<dof_id_type, std::set<Elem *, CompareDofObjectsByID>>
450 
452  std::set<SubdomainID> _secondary_ip_sub_ids;
453 
455  std::set<SubdomainID> _primary_ip_sub_ids;
456 
462  void projectSecondaryNodesSinglePair(SubdomainID lower_dimensional_primary_subdomain_id,
463  SubdomainID lower_dimensional_secondary_subdomain_id);
464 
468  void projectPrimaryNodesSinglePair(SubdomainID lower_dimensional_primary_subdomain_id,
469  SubdomainID lower_dimensional_secondary_subdomain_id);
470 
474  void
475  householderOrthogolization(const Point & normal, Point & tangent_one, Point & tangent_two) const;
476 
478  const bool _debug;
479 
481  const bool _on_displaced;
482 
484  const bool _periodic;
485 
487  const bool _distributed;
488 
491 
496 
500 
506 
508  std::unique_ptr<InputParameters> _output_params;
509 
512 };
513 
514 inline const std::pair<BoundaryID, BoundaryID> &
516 {
517  mooseAssert(_primary_secondary_boundary_id_pairs.size() == 1,
518  "We currently only support a single boundary pair per mortar generation object");
519 
521 }
std::set< SubdomainID > _primary_boundary_subdomain_ids
const std::pair< BoundaryID, BoundaryID > & primarySecondaryBoundaryIDPair() const
const Elem * getSecondaryLowerdElemFromSecondaryElem(dof_id_type secondary_elem_id) const
Return lower dimensional secondary element given its interior parent.
std::unordered_set< const Elem * > _inactive_local_lm_elems
List of inactive lagrange multiplier nodes (for elemental variables)
std::unordered_map< dof_id_type, std::set< Elem *, CompareDofObjectsByID > >::const_iterator MortarFilterIter
void computeInactiveLMNodes()
Get list of secondary nodes that don&#39;t contribute to interaction with any primary element...
std::vector< std::pair< BoundaryID, BoundaryID > > _primary_secondary_boundary_id_pairs
A list of primary/secondary boundary id pairs corresponding to each side of the mortar interface...
void clear()
Clears the mortar segment mesh and accompanying data structures.
std::set< SubdomainID > _secondary_boundary_subdomain_ids
The secondary/primary lower-dimensional boundary subdomain ids are the secondary/primary boundary ids...
std::unordered_map< dof_id_type, std::set< Elem *, CompareDofObjectsByID > > _secondary_elems_to_mortar_segments
We maintain a mapping from lower-dimensional secondary elements in the original mesh to (sets of) ele...
std::map< unsigned int, unsigned int > getPrimaryIpToLowerElementMap(const Elem &primary_elem, const Elem &primary_elem_ip, const Elem &lower_secondary_elem) const
Compute on-the-fly mapping from primary interior parent nodes to its corresponding lower dimensional ...
const bool _periodic
Whether this object will be generating a mortar segment mesh for periodic constraints.
Real _newton_tolerance
Newton solve tolerance for node projections.
TestClass subdomain_id_type
std::unordered_map< const Elem *, MortarSegmentInfo > _msm_elem_to_info
Map between Elems in the mortar segment mesh and their info structs.
const bool _distributed
Whether the mortar segment mesh is distributed.
Base class for MOOSE-based applications.
Definition: MooseApp.h:96
void buildMortarSegmentMesh()
Builds the mortar segment mesh once the secondary and primary node projections have been completed...
const std::set< SubdomainID > & primaryIPSubIDs() const
void msmStatistics()
Outputs mesh statistics for mortar segment mesh.
const std::unordered_set< const Elem * > & getInactiveLMElems() const
void buildCouplingInformation()
build the _mortar_interface_coupling data
std::vector< Point > getNodalNormals(const Elem &secondary_elem) const
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
static const std::string system_name
The name of the nodal normals system.
void projectSecondaryNodes()
Project secondary nodes (find xi^(2) values) to the closest points on the primary surface...
const std::unordered_map< dof_id_type, std::set< Elem *, CompareDofObjectsByID > > & secondariesToMortarSegments() const
const std::set< SubdomainID > & secondaryIPSubIDs() const
This class is a container/interface for the objects involved in automatic generation of mortar spaces...
const bool _debug
Whether to print debug output.
std::set< SubdomainID > _secondary_ip_sub_ids
All the secondary interior parent subdomain IDs associated with the mortar mesh.
void projectPrimaryNodes()
(Inverse) project primary nodes to the points on the secondary surface where they would have come fro...
int8_t boundary_id_type
AutomaticMortarGeneration(MooseApp &app, MeshBase &mesh_in, const std::pair< BoundaryID, BoundaryID > &boundary_key, const std::pair< SubdomainID, SubdomainID > &subdomain_key, bool on_displaced, bool periodic, const bool debug, const bool correct_edge_dropping, const Real minimum_projection_angle)
Must be constructed with a reference to the Mesh we are generating mortar spaces for.
void buildNodeToElemMaps()
Once the secondary_requested_boundary_ids and primary_requested_boundary_ids containers have been fil...
std::unordered_map< dof_id_type, const Elem * > _secondary_element_to_secondary_lowerd_element
Map from full dimensional secondary element id to lower dimensional secondary element.
std::map< std::tuple< dof_id_type, const Node *, const Elem * >, std::pair< Real, const Elem * > > _primary_node_and_elem_to_xi1_secondary_elem
Same type of container, but for mapping (Primary Node ID, Primary Node, Primary Elem) -> (xi^(1)...
boundary_id_type BoundaryID
void projectPrimaryNodesSinglePair(SubdomainID lower_dimensional_primary_subdomain_id, SubdomainID lower_dimensional_secondary_subdomain_id)
Helper function used internally by AutomaticMortarGeneration::project_primary_nodes().
std::unordered_map< dof_id_type, std::vector< const Elem * > > _nodes_to_primary_elem_map
const bool _correct_edge_dropping
Flag to enable regressed treatment of edge dropping where all LM DoFs on edge dropping element are st...
An inteface for the _console for outputting to the Console object.
std::array< MooseUtils::SemidynamicVector< Point, 9 >, 2 > getNodalTangents(const Elem &secondary_elem) const
Compute the two nodal tangents, which are built on-the-fly.
std::unique_ptr< InputParameters > _output_params
Storage for the input parameters used by the mortar nodal geometry output.
const std::unordered_map< dof_id_type, std::unordered_set< dof_id_type > > & mortarInterfaceCoupling() const
std::set< BoundaryID > _primary_requested_boundary_ids
The boundary ids corresponding to all the primary surfaces.
std::map< unsigned int, unsigned int > getSecondaryIpToLowerElementMap(const Elem &lower_secondary_elem) const
Compute on-the-fly mapping from secondary interior parent nodes to lower dimensional nodes...
subdomain_id_type SubdomainID
void computeNodalGeometry()
Computes and stores the nodal normal/tangent vectors in a local data structure instead of using the E...
Real _xi_tolerance
Tolerance for checking projection xi values.
void computeIncorrectEdgeDroppingInactiveLMNodes()
Computes inactive secondary nodes when incorrect edge dropping behavior is enabled (any node touching...
void buildMortarSegmentMesh3d()
Builds the mortar segment mesh once the secondary and primary node projections have been completed...
const std::unordered_map< dof_id_type, std::vector< const Elem * > > & nodesToSecondaryElem() const
bool onDisplaced() const
returns whether this object is on the displaced mesh
void projectSecondaryNodesSinglePair(SubdomainID lower_dimensional_primary_subdomain_id, SubdomainID lower_dimensional_secondary_subdomain_id)
Helper function responsible for projecting secondary nodes onto primary elements for a single primary...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unordered_map< std::pair< const Node *, const Elem * >, std::pair< Real, const Elem * > > _secondary_node_and_elem_to_xi2_primary_elem
Similar to the map above, but associates a (Secondary Node, Secondary Elem) pair to a (xi^(2)...
const bool _on_displaced
Whether this object is on the displaced mesh.
void initOutput()
initialize mortar-mesh based output
const std::unordered_map< const Elem *, MortarSegmentInfo > & mortarSegmentMeshElemToInfo() const
std::unordered_map< const Node *, std::array< Point, 2 > > _secondary_node_to_hh_nodal_tangents
Container for storing the nodal tangent/binormal vectors associated with each secondary node (Househo...
unsigned int mesh_dimension() const
const Real _minimum_projection_angle
Parameter to control which angle (in degrees) is admissible for the creation of mortar segments...
MeshBase & _mesh
Reference to the mesh stored in equation_systems.
const MeshBase & mortarSegmentMesh() const
std::vector< std::pair< SubdomainID, SubdomainID > > _primary_secondary_subdomain_id_pairs
A list of primary/secondary subdomain id pairs corresponding to each side of the mortar interface...
const std::unordered_set< const Node * > & getInactiveLMNodes() const
std::set< BoundaryID > _secondary_requested_boundary_ids
The boundary ids corresponding to all the secondary surfaces.
std::unordered_map< const Node *, Point > _secondary_node_to_nodal_normal
Container for storing the nodal normal vector associated with each secondary node.
std::unordered_set< const Node * > _inactive_local_lm_nodes
std::unordered_map< dof_id_type, std::unordered_set< dof_id_type > > _mortar_interface_coupling
Used by the AugmentSparsityOnInterface functor to determine whether a given Elem is coupled to any ot...
std::set< SubdomainID > _primary_ip_sub_ids
All the primary interior parent subdomain IDs associated with the mortar mesh.
void householderOrthogolization(const Point &normal, Point &tangent_one, Point &tangent_two) const
Householder orthogonalization procedure to obtain proper basis for tangent and binormal vectors...
std::vector< Point > getNormals(const Elem &secondary_elem, const std::vector< Point > &xi1_pts) const
Compute the normals at given reference points on a secondary element.
MooseApp & _app
The Moose app.
std::unique_ptr< MeshBase > _mortar_segment_mesh
1D Mesh of mortar segment elements which gets built by the call to build_mortar_segment_mesh().
uint8_t dof_id_type
std::unordered_map< dof_id_type, std::vector< const Elem * > > _nodes_to_secondary_elem_map
Map from nodes to connected lower-dimensional elements on the secondary/primary subdomains.
std::unordered_map< const Elem *, unsigned int > _lower_elem_to_side_id
Keeps track of the mapping between lower-dimensional elements and the side_id of the interior_parent ...
void computeInactiveLMElems()
Get list of secondary elems without any corresponding primary elements.