https://mooseframework.inl.gov
MortarConsumerInterface.C
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 // MOOSE Includes
12 #include "InputParameters.h"
13 #include "MooseObject.h"
14 #include "FEProblemBase.h"
15 #include "MooseMesh.h"
17 #include "Assembly.h"
19 #include "libmesh/quadrature.h"
20 
21 #include <algorithm>
22 
25 {
26  // Create InputParameters object that will be appended to the parameters for the inheriting object
28  // On a displaced mesh this will geometrically and algebraically ghost the entire interface
30  "AugmentSparsityOnInterface",
32  [](const InputParameters & obj_params, InputParameters & rm_params)
33  {
34  rm_params.set<bool>("use_displaced_mesh") = obj_params.get<bool>("use_displaced_mesh");
35  rm_params.set<BoundaryName>("secondary_boundary") =
36  obj_params.get<BoundaryName>("secondary_boundary");
37  rm_params.set<BoundaryName>("primary_boundary") =
38  obj_params.get<BoundaryName>("primary_boundary");
39  rm_params.set<SubdomainName>("secondary_subdomain") =
40  obj_params.get<SubdomainName>("secondary_subdomain");
41  rm_params.set<SubdomainName>("primary_subdomain") =
42  obj_params.get<SubdomainName>("primary_subdomain");
43  rm_params.set<bool>("ghost_point_neighbors") =
44  obj_params.get<bool>("ghost_point_neighbors");
45  rm_params.set<bool>("ghost_higher_d_neighbors") =
46  obj_params.get<bool>("ghost_higher_d_neighbors");
47  });
48 
49  params.addRequiredParam<BoundaryName>("primary_boundary",
50  "The name of the primary boundary sideset.");
51  params.addRequiredParam<BoundaryName>("secondary_boundary",
52  "The name of the secondary boundary sideset.");
53  params.addRequiredParam<SubdomainName>("primary_subdomain", "The name of the primary subdomain.");
54  params.addRequiredParam<SubdomainName>("secondary_subdomain",
55  "The name of the secondary subdomain.");
56  params.addParam<bool>(
57  "periodic",
58  false,
59  "Whether this constraint is going to be used to enforce a periodic condition. This has the "
60  "effect of changing the normals vector for projection from outward to inward facing");
61 
62  params.addParam<bool>(
63  "debug_mesh",
64  false,
65  "Whether this constraint is going to enable mortar segment mesh debug information. An exodus"
66  "file will be generated if the user sets this flag to true");
67 
68  params.addParam<bool>(
69  "correct_edge_dropping",
70  false,
71  "Whether to enable correct edge dropping treatment for mortar constraints. When disabled "
72  "any Lagrange Multiplier degree of freedom on a secondary element without full primary "
73  "contributions will be set (strongly) to 0.");
74 
75  params.addParam<bool>(
76  "interpolate_normals",
77  true,
78  "Whether to interpolate the nodal normals (e.g. classic idea of evaluating field at "
79  "quadrature points). If this is set to false, then non-interpolated nodal normals will be "
80  "used, and then the _normals member should be indexed with _i instead of _qp");
81 
82  params.addParam<bool>("ghost_point_neighbors",
83  false,
84  "Whether we should ghost point neighbors of secondary face elements, and "
85  "consequently also their mortar interface couples.");
86  params.addParam<Real>(
87  "minimum_projection_angle",
88  40.0,
89  "Parameter to control which angle (in degrees) is admissible for the creation of mortar "
90  "segments. If set to a value close to zero, very oblique projections are allowed, which "
91  "can result in mortar segments solving physics not meaningfully, and overprojection of "
92  "primary nodes onto the mortar segment mesh in extreme cases. This parameter is mostly "
93  "intended for mortar mesh debugging purposes in two dimensions.");
94 
95  params.addParam<bool>(
96  "ghost_higher_d_neighbors",
97  false,
98  "Whether we should ghost higher-dimensional neighbors. This is necessary when we are doing "
99  "second order mortar with finite volume primal variables, because in order for the method to "
100  "be second order we must use cell gradients, which couples in the neighbor cells.");
101 
102  return params;
103 }
104 
105 // Standard constructor
107  : _mci_fe_problem(*moose_object->getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
108  _mci_subproblem(*moose_object->getCheckedPointerParam<SubProblem *>("_subproblem")),
109  _mci_tid(moose_object->getParam<THREAD_ID>("_tid")),
110  _mci_mesh(_mci_subproblem.mesh()),
111  // all geometric assembly information should be correct for nl system number 0
112  _mci_assembly(_mci_subproblem.assembly(_mci_tid, 0)),
113  _mortar_data(_mci_fe_problem.mortarData()),
114  _secondary_id(
115  _mci_mesh.getBoundaryID(moose_object->getParam<BoundaryName>("secondary_boundary"))),
116  _primary_id(_mci_mesh.getBoundaryID(moose_object->getParam<BoundaryName>("primary_boundary"))),
117  _secondary_subdomain_id(
118  _mci_mesh.getSubdomainID(moose_object->getParam<SubdomainName>("secondary_subdomain"))),
119  _primary_subdomain_id(
120  _mci_mesh.getSubdomainID(moose_object->getParam<SubdomainName>("primary_subdomain"))),
121  _secondary_set({_secondary_id}),
122  _interpolate_normals(moose_object->getParam<bool>("interpolate_normals")),
123  _phys_points_secondary(_mci_assembly.qPointsFace()),
124  _phys_points_primary(_mci_assembly.qPointsFaceNeighbor()),
125  _qrule_msm(_mci_assembly.qRuleMortar()),
126  _qrule_face(_mci_assembly.qRuleFace()),
127  _lower_secondary_elem(_mci_assembly.lowerDElem()),
128  _lower_primary_elem(_mci_assembly.neighborLowerDElem()),
129  _JxW_msm(_mci_assembly.jxWMortar()),
130  _msm_elem(_mci_assembly.msmElem())
131 {
132  const bool displaced = moose_object->isParamValid("use_displaced_mesh")
133  ? moose_object->getParam<bool>("use_displaced_mesh")
134  : false;
135 
136  // Create the mortar interface if it hasn't already been created
137  _mci_fe_problem.createMortarInterface(
138  std::make_pair(_primary_id, _secondary_id),
139  std::make_pair(_primary_subdomain_id, _secondary_subdomain_id),
140  displaced,
141  moose_object->getParam<bool>("periodic"),
142  moose_object->getParam<bool>("debug_mesh"),
143  moose_object->getParam<bool>("correct_edge_dropping"),
144  moose_object->getParam<Real>("minimum_projection_angle"));
145 
146  _amg = &_mci_fe_problem.getMortarInterface(
147  std::make_pair(_primary_id, _secondary_id),
148  std::make_pair(_primary_subdomain_id, _secondary_subdomain_id),
149  displaced);
150 
151  const auto & secondary_set = _mortar_data.getHigherDimSubdomainIDs(_secondary_subdomain_id);
152  const auto & primary_set = _mortar_data.getHigherDimSubdomainIDs(_primary_subdomain_id);
153 
154  std::set_union(secondary_set.begin(),
155  secondary_set.end(),
156  primary_set.begin(),
157  primary_set.end(),
158  std::inserter(_higher_dim_subdomain_ids, _higher_dim_subdomain_ids.begin()));
159  _boundary_ids = {_secondary_id, _primary_id};
160 }
161 
162 void
164 {
165  if (interpolateNormals())
167  else
169 }
170 
171 void
173  ADReal & dual_number)
174 {
175  auto md_it = dual_number.derivatives().nude_data().begin();
176  auto mi_it = dual_number.derivatives().nude_indices().begin();
177 
178  auto d_it = dual_number.derivatives().nude_data().begin();
179 
180  for (auto i_it = dual_number.derivatives().nude_indices().begin();
181  i_it != dual_number.derivatives().nude_indices().end();
182  ++i_it, ++d_it)
183  if (*i_it != remove_derivative_index)
184  {
185  *mi_it = *i_it;
186  *md_it = *d_it;
187  ++mi_it;
188  ++md_it;
189  }
190 
191  std::size_t n_indices = md_it - dual_number.derivatives().nude_data().begin();
192  dual_number.derivatives().nude_indices().resize(n_indices);
193  dual_number.derivatives().nude_data().resize(n_indices);
194 }
const BoundaryID _secondary_id
Boundary ID for the secondary surface.
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
MortarConsumerInterface(const MooseObject *moose_object)
const libMesh::QBase *const & _qrule_face
The arbitrary quadrature rule on the lower dimensional secondary face.
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
MeshBase & mesh
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
std::vector< Point > getNodalNormals(const Elem &secondary_elem) const
void addRelationshipManager(const std::string &name, Moose::RelationshipManagerType rm_type, Moose::RelationshipManagerInputParameterCallback input_parameter_callback=nullptr)
Tells MOOSE about a RelationshipManager that this object needs.
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:42
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
InputParameters emptyInputParameters()
SubdomainID getSubdomainID(const SubdomainName &subdomain_name, const MeshBase &mesh)
Gets the subdomain ID associated with the given SubdomainName.
BoundaryID getBoundaryID(const BoundaryName &boundary_name, const MeshBase &mesh)
Gets the boundary ID associated with the given BoundaryName.
Elem const *const & _lower_secondary_elem
The secondary face lower dimensional element (not the mortar element!).
const AutomaticMortarGeneration & amg() const
Retrieve the automatic mortar generation object associated with this constraint.
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
std::vector< Point > _normals
the normals
bool interpolateNormals() const
Whether to interpolate the nodal normals (e.g.
static void trimDerivative(dof_id_type remove_derivative_index, ADReal &dual_number)
Get rid of AD derivative entries by dof index.
void setNormals()
Set the normals vector.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:78
const std::vector< Point > & get_points() const
static InputParameters validParams()
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
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.
unsigned int THREAD_ID
Definition: MooseTypes.h:237
uint8_t dof_id_type