This class is a container/interface for the objects involved in automatic generation of mortar spaces. More...
#include <AutomaticMortarGeneration.h>
Public Types | |
using | MortarFilterIter = std::unordered_map< dof_id_type, std::set< Elem *, CompareDofObjectsByID > >::const_iterator |
Public Member Functions | |
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. More... | |
void | buildNodeToElemMaps () |
Once the secondary_requested_boundary_ids and primary_requested_boundary_ids containers have been filled in, call this function to build node-to-Elem maps for the lower-dimensional elements. More... | |
void | computeNodalGeometry () |
Computes and stores the nodal normal/tangent vectors in a local data structure instead of using the ExplicitSystem/NumericVector approach. More... | |
void | projectSecondaryNodes () |
Project secondary nodes (find xi^(2) values) to the closest points on the primary surface. More... | |
void | projectPrimaryNodes () |
(Inverse) project primary nodes to the points on the secondary surface where they would have come from (find (xi^(1) values)). More... | |
void | buildMortarSegmentMesh () |
Builds the mortar segment mesh once the secondary and primary node projections have been completed. More... | |
void | buildMortarSegmentMesh3d () |
Builds the mortar segment mesh once the secondary and primary node projections have been completed. More... | |
void | msmStatistics () |
Outputs mesh statistics for mortar segment mesh. More... | |
void | clear () |
Clears the mortar segment mesh and accompanying data structures. More... | |
bool | onDisplaced () const |
returns whether this object is on the displaced mesh More... | |
std::vector< Point > | getNodalNormals (const Elem &secondary_elem) const |
std::array< MooseUtils::SemidynamicVector< Point, 9 >, 2 > | getNodalTangents (const Elem &secondary_elem) const |
Compute the two nodal tangents, which are built on-the-fly. More... | |
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. More... | |
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 nodes. More... | |
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. More... | |
std::vector< Point > | getNormals (const Elem &secondary_elem, const std::vector< Real > &oned_xi1_pts) const |
Compute the normals at given reference points on a secondary element. More... | |
const Elem * | getSecondaryLowerdElemFromSecondaryElem (dof_id_type secondary_elem_id) const |
Return lower dimensional secondary element given its interior parent. More... | |
void | computeInactiveLMNodes () |
Get list of secondary nodes that don't contribute to interaction with any primary element. More... | |
void | computeIncorrectEdgeDroppingInactiveLMNodes () |
Computes inactive secondary nodes when incorrect edge dropping behavior is enabled (any node touching a partially or fully dropped element is dropped) More... | |
void | computeInactiveLMElems () |
Get list of secondary elems without any corresponding primary elements. More... | |
const std::unordered_map< dof_id_type, std::unordered_set< dof_id_type > > & | mortarInterfaceCoupling () const |
const std::pair< BoundaryID, BoundaryID > & | primarySecondaryBoundaryIDPair () const |
const MeshBase & | mortarSegmentMesh () const |
const std::unordered_map< const Elem *, MortarSegmentInfo > & | mortarSegmentMeshElemToInfo () const |
int | dim () const |
const std::unordered_set< const Node * > & | getInactiveLMNodes () const |
const std::unordered_set< const Elem * > & | getInactiveLMElems () const |
bool | incorrectEdgeDropping () const |
std::vector< MortarFilterIter > | secondariesToMortarSegments (const Node &node) const |
const std::unordered_map< dof_id_type, std::set< Elem *, CompareDofObjectsByID > > & | secondariesToMortarSegments () const |
const std::set< SubdomainID > & | secondaryIPSubIDs () const |
const std::set< SubdomainID > & | primaryIPSubIDs () const |
const std::unordered_map< dof_id_type, std::vector< const Elem * > > & | nodesToSecondaryElem () const |
void | initOutput () |
initialize mortar-mesh based output More... | |
Public Attributes | |
const ConsoleStream | _console |
An instance of helper class to write streams to the Console objects. More... | |
Static Public Attributes | |
static const std::string | system_name |
The name of the nodal normals system. More... | |
Private Member Functions | |
void | buildCouplingInformation () |
build the _mortar_interface_coupling data More... | |
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/secondary pair. More... | |
void | projectPrimaryNodesSinglePair (SubdomainID lower_dimensional_primary_subdomain_id, SubdomainID lower_dimensional_secondary_subdomain_id) |
Helper function used internally by AutomaticMortarGeneration::project_primary_nodes(). More... | |
void | householderOrthogolization (const Point &normal, Point &tangent_one, Point &tangent_two) const |
Householder orthogonalization procedure to obtain proper basis for tangent and binormal vectors. More... | |
Private Attributes | |
MooseApp & | _app |
The Moose app. More... | |
MeshBase & | _mesh |
Reference to the mesh stored in equation_systems. More... | |
std::set< BoundaryID > | _secondary_requested_boundary_ids |
The boundary ids corresponding to all the secondary surfaces. More... | |
std::set< BoundaryID > | _primary_requested_boundary_ids |
The boundary ids corresponding to all the primary surfaces. More... | |
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. More... | |
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. More... | |
std::unordered_map< dof_id_type, std::vector< const Elem * > > | _nodes_to_primary_elem_map |
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), primary Elem) pair. More... | |
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), Secondary Elem) where they are inverse-projected along the nodal normal direction. More... | |
std::unique_ptr< MeshBase > | _mortar_segment_mesh |
1D Mesh of mortar segment elements which gets built by the call to build_mortar_segment_mesh(). More... | |
std::unordered_map< const Elem *, MortarSegmentInfo > | _msm_elem_to_info |
Map between Elems in the mortar segment mesh and their info structs. More... | |
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 which they are. More... | |
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. More... | |
std::set< SubdomainID > | _secondary_boundary_subdomain_ids |
The secondary/primary lower-dimensional boundary subdomain ids are the secondary/primary boundary ids. More... | |
std::set< SubdomainID > | _primary_boundary_subdomain_ids |
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 others across the gap, and to explicitly set up the dependence between interior_parent() elements on the secondary side and their lower-dimensional sides which are on the interface. More... | |
std::unordered_map< const Node *, Point > | _secondary_node_to_nodal_normal |
Container for storing the nodal normal vector associated with each secondary node. More... | |
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 (Householder approach). More... | |
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. More... | |
std::unordered_set< const Node * > | _inactive_local_lm_nodes |
std::unordered_set< const Elem * > | _inactive_local_lm_elems |
List of inactive lagrange multiplier nodes (for elemental variables) More... | |
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) elements in mortar_segment_mesh. More... | |
std::set< SubdomainID > | _secondary_ip_sub_ids |
All the secondary interior parent subdomain IDs associated with the mortar mesh. More... | |
std::set< SubdomainID > | _primary_ip_sub_ids |
All the primary interior parent subdomain IDs associated with the mortar mesh. More... | |
const bool | _debug |
Whether to print debug output. More... | |
const bool | _on_displaced |
Whether this object is on the displaced mesh. More... | |
const bool | _periodic |
Whether this object will be generating a mortar segment mesh for periodic constraints. More... | |
const bool | _distributed |
Whether the mortar segment mesh is distributed. More... | |
Real | _newton_tolerance = 1e-12 |
Newton solve tolerance for node projections. More... | |
Real | _xi_tolerance = 1e-6 |
Tolerance for checking projection xi values. More... | |
const bool | _correct_edge_dropping |
Flag to enable regressed treatment of edge dropping where all LM DoFs on edge dropping element are strongly set to 0. More... | |
const Real | _minimum_projection_angle |
Parameter to control which angle (in degrees) is admissible for the creation of mortar segments. More... | |
std::unique_ptr< InputParameters > | _output_params |
Storage for the input parameters used by the mortar nodal geometry output. More... | |
Friends | |
class | MortarNodalGeometryOutput |
class | AugmentSparsityOnInterface |
This class is a container/interface for the objects involved in automatic generation of mortar spaces.
Definition at line 56 of file AutomaticMortarGeneration.h.
using AutomaticMortarGeneration::MortarFilterIter = std::unordered_map<dof_id_type, std::set<Elem *, CompareDofObjectsByID> >::const_iterator |
Definition at line 300 of file AutomaticMortarGeneration.h.
AutomaticMortarGeneration::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.
Definition at line 213 of file AutomaticMortarGeneration.C.
|
private |
build the _mortar_interface_coupling
data
Definition at line 1351 of file AutomaticMortarGeneration.C.
Referenced by buildMortarSegmentMesh(), and buildMortarSegmentMesh3d().
void AutomaticMortarGeneration::buildMortarSegmentMesh | ( | ) |
Builds the mortar segment mesh once the secondary and primary node projections have been completed.
Inputs:
Outputs:
Defined in the file build_mortar_segment_mesh.C.
This was a change to how inactive LM DoFs are handled. Now mortar segment elements are not used in assembly if there is no corresponding primary element and inactive LM DoFs (those with no contribution to an active primary element) are zeroed.
Definition at line 434 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
void AutomaticMortarGeneration::buildMortarSegmentMesh3d | ( | ) |
Builds the mortar segment mesh once the secondary and primary node projections have been completed.
Inputs:
Outputs:
Step 1: Build mortar segments for all secondary elements
Step 1.1: Linearize secondary face elements
For first order face elements (Tri3 and Quad4) elements are simply linearized around center For second order (Tri6 and Quad9) and third order (Tri7) face elements, elements are sub-divided into four first order elements then each of the sub-elements is linearized around their respective centers For Quad8 elements, they are sub-divided into one quad and four triangle elements and each sub-element is linearized around their respective centers
Step 1.2: Coarse screening using a k-d tree to find nodes on the primary interface that are 'close to' a center point of the secondary element.
Step 1.3: Loop through primary candidate nodes, create mortar segments
Once an element with non-trivial projection onto secondary element identified, switch to breadth-first search (drop all current candidates and add only neighbors of elements with non-trivial overlap)
Step 1.3.2: Sub-divide primary element candidate, then project onto secondary sub-elements, perform polygon clipping, and triangulate to form mortar segments
Step 1.3.3: Create mortar segments and add to mortar segment mesh
Definition at line 931 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
void AutomaticMortarGeneration::buildNodeToElemMaps | ( | ) |
Once the secondary_requested_boundary_ids and primary_requested_boundary_ids containers have been filled in, call this function to build node-to-Elem maps for the lower-dimensional elements.
Definition at line 285 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
void AutomaticMortarGeneration::clear | ( | ) |
Clears the mortar segment mesh and accompanying data structures.
Definition at line 266 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
void AutomaticMortarGeneration::computeInactiveLMElems | ( | ) |
Get list of secondary elems without any corresponding primary elements.
Used to enforce zero values on inactive DoFs of elemental variables.
Definition at line 1662 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
void AutomaticMortarGeneration::computeInactiveLMNodes | ( | ) |
Get list of secondary nodes that don't contribute to interaction with any primary element.
Used to enforce zero values on inactive DoFs of nodal variables.
Definition at line 1578 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
void AutomaticMortarGeneration::computeIncorrectEdgeDroppingInactiveLMNodes | ( | ) |
Computes inactive secondary nodes when incorrect edge dropping behavior is enabled (any node touching a partially or fully dropped element is dropped)
Definition at line 1486 of file AutomaticMortarGeneration.C.
Referenced by computeInactiveLMNodes().
void AutomaticMortarGeneration::computeNodalGeometry | ( | ) |
Computes and stores the nodal normal/tangent vectors in a local data structure instead of using the ExplicitSystem/NumericVector approach.
This design was triggered by the way that the GhostingFunctor operates, but I think it is a better/more efficient way to do it anyway.
The _periodic flag tells us whether we want to inward vs outward facing normals
Definition at line 1710 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
|
inline |
Definition at line 279 of file AutomaticMortarGeneration.h.
Referenced by computeInactiveLMElems(), computeIncorrectEdgeDroppingInactiveLMNodes(), computeNodalGeometry(), Moose::Mortar::loopOverMortarSegments(), and MortarData::update().
|
inline |
Definition at line 292 of file AutomaticMortarGeneration.h.
Referenced by ComputeMortarFunctor::operator()().
|
inline |
Definition at line 284 of file AutomaticMortarGeneration.h.
Referenced by ComputeMortarFunctor::operator()().
std::vector< Point > AutomaticMortarGeneration::getNodalNormals | ( | const Elem & | secondary_elem | ) | const |
secondary_elem
Definition at line 323 of file AutomaticMortarGeneration.C.
Referenced by buildMortarSegmentMesh3d(), getNormals(), and MortarConsumerInterface::setNormals().
std::array< MooseUtils::SemidynamicVector< Point, 9 >, 2 > AutomaticMortarGeneration::getNodalTangents | ( | const Elem & | secondary_elem | ) | const |
Compute the two nodal tangents, which are built on-the-fly.
secondary_elem
Definition at line 376 of file AutomaticMortarGeneration.C.
std::vector< Point > AutomaticMortarGeneration::getNormals | ( | const Elem & | secondary_elem, |
const std::vector< Point > & | xi1_pts | ||
) | const |
Compute the normals at given reference points on a secondary element.
secondary_elem | The secondary element used to query for associated nodal normals |
xi1_pts | The reference points on the secondary element to evaluate the normals at. The points should only be non-zero in the zeroth entry because right now our mortar mesh elements are always 1D |
Definition at line 405 of file AutomaticMortarGeneration.C.
Referenced by buildMortarSegmentMesh(), getNormals(), and MortarConsumerInterface::setNormals().
std::vector< Point > AutomaticMortarGeneration::getNormals | ( | const Elem & | secondary_elem, |
const std::vector< Real > & | oned_xi1_pts | ||
) | const |
Compute the normals at given reference points on a secondary element.
secondary_elem | The secondary element used to query for associated nodal normals |
1d_xi1_pts | The reference points on the secondary element to evaluate the normals at. The "points" are single reals corresponding to xi because right now our mortar mesh elements are always 1D |
Definition at line 394 of file AutomaticMortarGeneration.C.
std::map< unsigned int, unsigned int > AutomaticMortarGeneration::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 nodes.
Definition at line 359 of file AutomaticMortarGeneration.C.
std::map< unsigned int, unsigned int > AutomaticMortarGeneration::getSecondaryIpToLowerElementMap | ( | const Elem & | lower_secondary_elem | ) | const |
Compute on-the-fly mapping from secondary interior parent nodes to lower dimensional nodes.
Definition at line 343 of file AutomaticMortarGeneration.C.
const Elem * AutomaticMortarGeneration::getSecondaryLowerdElemFromSecondaryElem | ( | dof_id_type | secondary_elem_id | ) | const |
Return lower dimensional secondary element given its interior parent.
Helpful outside the mortar generation to locate mortar-related quantities.
secondary_elem_id | The secondary interior parent element id used to query for associated lower dimensional element |
Definition at line 333 of file AutomaticMortarGeneration.C.
|
private |
Householder orthogonalization procedure to obtain proper basis for tangent and binormal vectors.
Definition at line 1803 of file AutomaticMortarGeneration.C.
Referenced by computeNodalGeometry().
|
inline |
Definition at line 297 of file AutomaticMortarGeneration.h.
Referenced by ComputeMortarFunctor::operator()().
void AutomaticMortarGeneration::initOutput | ( | ) |
initialize mortar-mesh based output
Definition at line 247 of file AutomaticMortarGeneration.C.
|
inline |
Definition at line 256 of file AutomaticMortarGeneration.h.
|
inline |
Definition at line 269 of file AutomaticMortarGeneration.h.
Referenced by Moose::Mortar::loopOverMortarSegments().
|
inline |
Definition at line 274 of file AutomaticMortarGeneration.h.
Referenced by Moose::Mortar::loopOverMortarSegments().
void AutomaticMortarGeneration::msmStatistics | ( | ) |
Outputs mesh statistics for mortar segment mesh.
Definition at line 1422 of file AutomaticMortarGeneration.C.
Referenced by buildMortarSegmentMesh3d().
|
inline |
Definition at line 333 of file AutomaticMortarGeneration.h.
|
inline |
returns whether this object is on the displaced mesh
Definition at line 171 of file AutomaticMortarGeneration.h.
|
inline |
Definition at line 328 of file AutomaticMortarGeneration.h.
Referenced by Moose::Mortar::setupMortarMaterials().
|
inline |
Definition at line 515 of file AutomaticMortarGeneration.h.
Referenced by Moose::Mortar::loopOverMortarSegments(), and Moose::Mortar::setupMortarMaterials().
void AutomaticMortarGeneration::projectPrimaryNodes | ( | ) |
(Inverse) project primary nodes to the points on the secondary surface where they would have come from (find (xi^(1) values)).
Inputs:
Outputs:
Defined in the file project_primary_nodes.C.
Definition at line 2165 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
|
private |
Helper function used internally by AutomaticMortarGeneration::project_primary_nodes().
Definition at line 2174 of file AutomaticMortarGeneration.C.
Referenced by projectPrimaryNodes().
void AutomaticMortarGeneration::projectSecondaryNodes | ( | ) |
Project secondary nodes (find xi^(2) values) to the closest points on the primary surface.
Inputs:
Outputs:
Defined in the file project_secondary_nodes.C.
Definition at line 1849 of file AutomaticMortarGeneration.C.
Referenced by MortarData::update().
|
private |
Helper function responsible for projecting secondary nodes onto primary elements for a single primary/secondary pair.
Called by the class member AutomaticMortarGeneration::project_secondary_nodes().
Definition at line 1858 of file AutomaticMortarGeneration.C.
Referenced by projectSecondaryNodes().
std::vector< AutomaticMortarGeneration::MortarFilterIter > AutomaticMortarGeneration::secondariesToMortarSegments | ( | const Node & | node | ) | const |
Definition at line 2378 of file AutomaticMortarGeneration.C.
Referenced by MortarUserObjectThread::operator()(), and ComputeMortarFunctor::operator()().
|
inline |
Definition at line 315 of file AutomaticMortarGeneration.h.
|
inline |
Definition at line 323 of file AutomaticMortarGeneration.h.
Referenced by Moose::Mortar::setupMortarMaterials().
|
friend |
Definition at line 511 of file AutomaticMortarGeneration.h.
|
friend |
Definition at line 510 of file AutomaticMortarGeneration.h.
|
private |
The Moose app.
Definition at line 350 of file AutomaticMortarGeneration.h.
Referenced by initOutput().
|
inherited |
An instance of helper class to write streams to the Console objects.
Definition at line 31 of file ConsoleStreamInterface.h.
Referenced by IterationAdaptiveDT::acceptStep(), MeshOnlyAction::act(), SetupDebugAction::act(), MaterialOutputAction::act(), Adaptivity::adaptMesh(), FEProblemBase::adaptMesh(), PerfGraph::addToExecutionList(), SimplePredictor::apply(), SystemBase::applyScalingFactors(), MultiApp::backup(), FEProblemBase::backupMultiApps(), CoarsenedPiecewiseLinear::buildCoarsenedGrid(), MeshDiagnosticsGenerator::checkElementOverlap(), MeshDiagnosticsGenerator::checkElementTypes(), MeshDiagnosticsGenerator::checkElementVolumes(), FEProblemBase::checkExceptionAndStopSolve(), SolverSystem::checkInvalidSolution(), MeshDiagnosticsGenerator::checkLocalJacobians(), MeshDiagnosticsGenerator::checkNonConformalMesh(), MeshDiagnosticsGenerator::checkNonConformalMeshFromAdaptivity(), MeshDiagnosticsGenerator::checkNonPlanarSides(), FEProblemBase::checkProblemIntegrity(), ReferenceResidualProblem::checkRelativeConvergence(), MeshDiagnosticsGenerator::checkSidesetsOrientation(), IterationAdaptiveDT::computeAdaptiveDT(), Transient::computeConstrainedDT(), FixedPointSolve::computeCustomConvergencePostprocessor(), NonlinearSystemBase::computeDamping(), IterationAdaptiveDT::computeDT(), IterationAdaptiveDT::computeFailedDT(), IterationAdaptiveDT::computeInitialDT(), IterationAdaptiveDT::computeInterpolationDT(), FEProblemBase::computeLinearSystemTags(), NonlinearSystemBase::computeScaling(), Problem::console(), IterationAdaptiveDT::constrainStep(), TimeStepper::constrainStep(), MultiApp::createApp(), FEProblemBase::execMultiApps(), FEProblemBase::execMultiAppTransfers(), MessageFromInput::execute(), Steady::execute(), Eigenvalue::execute(), ActionWarehouse::executeActionsWithAction(), ActionWarehouse::executeAllActions(), ElementQualityChecker::finalize(), FEProblemBase::finishMultiAppStep(), MeshRepairGenerator::fixOverlappingNodes(), CoarsenBlockGenerator::generate(), MeshGenerator::generateInternal(), VariableCondensationPreconditioner::getDofToCondense(), InversePowerMethod::init(), NonlinearEigen::init(), FEProblemBase::initialAdaptMesh(), EigenExecutionerBase::inversePowerIteration(), FEProblemBase::joinAndFinalize(), Transient::keepGoing(), IterationAdaptiveDT::limitDTByFunction(), IterationAdaptiveDT::limitDTToPostprocessorValue(), FEProblemBase::logAdd(), EigenExecutionerBase::makeBXConsistent(), Console::meshChanged(), MooseBaseErrorInterface::mooseDeprecated(), MooseBaseErrorInterface::mooseInfo(), MooseBaseErrorInterface::mooseWarning(), MooseBaseErrorInterface::mooseWarningNonPrefixed(), ReferenceResidualProblem::nonlinearConvergenceSetup(), ReporterDebugOutput::output(), PerfGraphOutput::output(), MaterialPropertyDebugOutput::output(), DOFMapOutput::output(), VariableResidualNormsDebugOutput::output(), Console::output(), ControlOutput::outputActiveObjects(), ControlOutput::outputChangedControls(), ControlOutput::outputControls(), Console::outputInput(), Console::outputPostprocessors(), PseudoTimestep::outputPseudoTimestep(), Console::outputReporters(), Console::outputScalarVariables(), Console::outputSystemInformation(), FEProblemBase::possiblyRebuildGeomSearchPatches(), EigenExecutionerBase::postExecute(), AB2PredictorCorrector::postSolve(), ActionWarehouse::printActionDependencySets(), SolutionInvalidity::printDebug(), EigenExecutionerBase::printEigenvalue(), SecantSolve::printFixedPointConvergenceHistory(), SteffensenSolve::printFixedPointConvergenceHistory(), PicardSolve::printFixedPointConvergenceHistory(), FixedPointSolve::printFixedPointConvergenceReason(), PerfGraphLivePrint::printLiveMessage(), MaterialPropertyDebugOutput::printMaterialMap(), PerfGraphLivePrint::printStats(), projectPrimaryNodesSinglePair(), projectSecondaryNodesSinglePair(), CoarsenBlockGenerator::recursiveCoarsen(), SolutionTimeAdaptiveDT::rejectStep(), MultiApp::restore(), FEProblemBase::restoreMultiApps(), NonlinearSystemBase::setInitialSolution(), Checkpoint::shouldOutput(), SubProblem::showFunctorRequestors(), SubProblem::showFunctors(), FullSolveMultiApp::showStatusMessage(), FEProblemSolve::solve(), FixedPointSolve::solve(), NonlinearSystem::solve(), EigenProblem::solve(), LStableDirk2::solve(), LStableDirk3::solve(), ImplicitMidpoint::solve(), ExplicitTVDRK2::solve(), LStableDirk4::solve(), AStableDirk4::solve(), ExplicitRK2::solve(), TransientMultiApp::solveStep(), FixedPointSolve::solveStep(), PerfGraphLivePrint::start(), AB2PredictorCorrector::step(), NonlinearEigen::takeStep(), Transient::takeStep(), Console::writeTimestepInformation(), Console::writeVariableNorms(), and FEProblemBase::~FEProblemBase().
|
private |
Flag to enable regressed treatment of edge dropping where all LM DoFs on edge dropping element are strongly set to 0.
Definition at line 499 of file AutomaticMortarGeneration.h.
Referenced by computeInactiveLMElems(), computeInactiveLMNodes(), and incorrectEdgeDropping().
|
private |
Whether to print debug output.
Definition at line 478 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), initOutput(), projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
Whether the mortar segment mesh is distributed.
Definition at line 487 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration().
|
private |
List of inactive lagrange multiplier nodes (for elemental variables)
Definition at line 443 of file AutomaticMortarGeneration.h.
Referenced by computeInactiveLMElems(), and getInactiveLMElems().
|
private |
Definition at line 440 of file AutomaticMortarGeneration.h.
Referenced by computeInactiveLMNodes(), computeIncorrectEdgeDroppingInactiveLMNodes(), and getInactiveLMNodes().
|
private |
Keeps track of the mapping between lower-dimensional elements and the side_id of the interior_parent which they are.
Definition at line 408 of file AutomaticMortarGeneration.h.
Referenced by clear().
|
private |
Reference to the mesh stored in equation_systems.
Definition at line 353 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), buildCouplingInformation(), buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), buildNodeToElemMaps(), computeInactiveLMElems(), computeInactiveLMNodes(), computeIncorrectEdgeDroppingInactiveLMNodes(), computeNodalGeometry(), dim(), getNormals(), msmStatistics(), projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
Parameter to control which angle (in degrees) is admissible for the creation of mortar segments.
If set to a value close to zero, very oblique projections are allowed, which can result in mortar segments solving physics not meaningfully and overprojection of primary nodes onto the mortar segment mesh in extreme cases. This parameter is mostly intended for mortar mesh debugging purposes in 2D.
Definition at line 505 of file AutomaticMortarGeneration.h.
Referenced by projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
Used by the AugmentSparsityOnInterface functor to determine whether a given Elem is coupled to any others across the gap, and to explicitly set up the dependence between interior_parent() elements on the secondary side and their lower-dimensional sides which are on the interface.
This latter type of coupling must be explicitly declared when there is no primary_elem for a given mortar segment and you are using e.g. a P^1-P^0 discretization which does not induce the coupling automatically.
Definition at line 427 of file AutomaticMortarGeneration.h.
Referenced by buildCouplingInformation(), clear(), and mortarInterfaceCoupling().
|
private |
1D Mesh of mortar segment elements which gets built by the call to build_mortar_segment_mesh().
Definition at line 399 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), clear(), computeInactiveLMElems(), computeInactiveLMNodes(), computeIncorrectEdgeDroppingInactiveLMNodes(), mortarSegmentMesh(), and msmStatistics().
|
private |
Map between Elems in the mortar segment mesh and their info structs.
This gets filled in by the call to build_mortar_segment_mesh().
Definition at line 404 of file AutomaticMortarGeneration.h.
Referenced by buildCouplingInformation(), buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), clear(), computeInactiveLMElems(), computeInactiveLMNodes(), computeIncorrectEdgeDroppingInactiveLMNodes(), and mortarSegmentMeshElemToInfo().
|
private |
Newton solve tolerance for node projections.
Definition at line 490 of file AutomaticMortarGeneration.h.
Referenced by projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
Definition at line 367 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), buildNodeToElemMaps(), clear(), projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
Map from nodes to connected lower-dimensional elements on the secondary/primary subdomains.
Definition at line 366 of file AutomaticMortarGeneration.h.
Referenced by buildNodeToElemMaps(), clear(), nodesToSecondaryElem(), projectPrimaryNodesSinglePair(), projectSecondaryNodesSinglePair(), and secondariesToMortarSegments().
|
private |
Whether this object is on the displaced mesh.
Definition at line 481 of file AutomaticMortarGeneration.h.
Referenced by initOutput(), and onDisplaced().
|
private |
Storage for the input parameters used by the mortar nodal geometry output.
Definition at line 508 of file AutomaticMortarGeneration.h.
Referenced by initOutput().
|
private |
Whether this object will be generating a mortar segment mesh for periodic constraints.
Definition at line 484 of file AutomaticMortarGeneration.h.
Referenced by computeNodalGeometry(), and getNormals().
|
private |
Definition at line 417 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), and buildNodeToElemMaps().
|
private |
All the primary interior parent subdomain IDs associated with the mortar mesh.
Definition at line 455 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), clear(), and primaryIPSubIDs().
|
private |
Same type of container, but for mapping (Primary Node ID, Primary Node, Primary Elem) -> (xi^(1), Secondary Elem) where they are inverse-projected along the nodal normal direction.
Note that the first item of the key, the primary node ID, is important for storing the key-value pairs in a consistent order across processes, e.g. this container has to be ordered!
Definition at line 395 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), clear(), projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
The boundary ids corresponding to all the primary surfaces.
Definition at line 359 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), and buildNodeToElemMaps().
|
private |
A list of primary/secondary boundary id pairs corresponding to each side of the mortar interface.
Definition at line 363 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), buildMortarSegmentMesh(), initOutput(), and primarySecondaryBoundaryIDPair().
|
private |
A list of primary/secondary subdomain id pairs corresponding to each side of the mortar interface.
Definition at line 412 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), buildMortarSegmentMesh3d(), computeInactiveLMElems(), computeInactiveLMNodes(), computeIncorrectEdgeDroppingInactiveLMNodes(), msmStatistics(), projectPrimaryNodes(), and projectSecondaryNodes().
|
private |
The secondary/primary lower-dimensional boundary subdomain ids are the secondary/primary boundary ids.
Definition at line 416 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), buildMortarSegmentMesh(), buildNodeToElemMaps(), and computeNodalGeometry().
|
private |
Map from full dimensional secondary element id to lower dimensional secondary element.
Definition at line 437 of file AutomaticMortarGeneration.h.
Referenced by clear(), computeNodalGeometry(), and getSecondaryLowerdElemFromSecondaryElem().
|
private |
We maintain a mapping from lower-dimensional secondary elements in the original mesh to (sets of) elements in mortar_segment_mesh.
This allows us to quickly determine which elements need to be split.
Definition at line 449 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), clear(), and secondariesToMortarSegments().
|
private |
All the secondary interior parent subdomain IDs associated with the mortar mesh.
Definition at line 452 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), buildMortarSegmentMesh3d(), clear(), and secondaryIPSubIDs().
|
private |
Similar to the map above, but associates a (Secondary Node, Secondary Elem) pair to a (xi^(2), primary Elem) pair.
This allows a single secondary node, which is potentially connected to two elements on the secondary side, to be associated with multiple primary Elem/xi^(2) values to handle the case where the primary and secondary nodes are "matching". In this configuration:
A B o--—o--—o (secondary orientation ->) | v ---—x---— (primary orientation <-) C D
The entries in the map should be: (Elem A, Node 1) -> (Elem C, xi^(2)=-1) (Elem B, Node 0) -> (Elem D, xi^(2)=+1)
Definition at line 387 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), clear(), and projectSecondaryNodesSinglePair().
|
private |
Container for storing the nodal tangent/binormal vectors associated with each secondary node (Householder approach).
Definition at line 434 of file AutomaticMortarGeneration.h.
Referenced by clear(), computeNodalGeometry(), and getNodalTangents().
|
private |
Container for storing the nodal normal vector associated with each secondary node.
Definition at line 430 of file AutomaticMortarGeneration.h.
Referenced by clear(), computeNodalGeometry(), getNodalNormals(), projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
private |
The boundary ids corresponding to all the secondary surfaces.
Definition at line 356 of file AutomaticMortarGeneration.h.
Referenced by AutomaticMortarGeneration(), and buildNodeToElemMaps().
|
private |
Tolerance for checking projection xi values.
Usually we are checking whether we projected onto a certain element (in which case -1 <= xi <= 1) or whether we should have already projected a primary node (in which case we error if abs(xi) is sufficiently close to 1)
Definition at line 495 of file AutomaticMortarGeneration.h.
Referenced by buildMortarSegmentMesh(), projectPrimaryNodesSinglePair(), and projectSecondaryNodesSinglePair().
|
static |
The name of the nodal normals system.
We store this in one place so it's easy to change later.
Definition at line 63 of file AutomaticMortarGeneration.h.