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"
16 #include "MortarData.h"
17 #include "Assembly.h"
18 #include "libmesh/quadrature.h"
19 
20 #include <algorithm>
21 
24 {
25  // Create InputParameters object that will be appended to the parameters for the inheriting object
27  // On a displaced mesh this will geometrically and algebraically ghost the entire interface
29  "AugmentSparsityOnInterface",
31  [](const InputParameters & obj_params, InputParameters & rm_params)
32  {
33  rm_params.set<bool>("use_displaced_mesh") = obj_params.get<bool>("use_displaced_mesh");
34  rm_params.set<BoundaryName>("secondary_boundary") =
35  obj_params.get<BoundaryName>("secondary_boundary");
36  rm_params.set<BoundaryName>("primary_boundary") =
37  obj_params.get<BoundaryName>("primary_boundary");
38  rm_params.set<SubdomainName>("secondary_subdomain") =
39  obj_params.get<SubdomainName>("secondary_subdomain");
40  rm_params.set<SubdomainName>("primary_subdomain") =
41  obj_params.get<SubdomainName>("primary_subdomain");
42  rm_params.set<bool>("ghost_point_neighbors") =
43  obj_params.get<bool>("ghost_point_neighbors");
44  rm_params.set<bool>("ghost_higher_d_neighbors") =
45  obj_params.get<bool>("ghost_higher_d_neighbors");
46  });
47 
48  params.addRequiredParam<BoundaryName>("primary_boundary",
49  "The name of the primary boundary sideset.");
50  params.addRequiredParam<BoundaryName>("secondary_boundary",
51  "The name of the secondary boundary sideset.");
52  params.addRequiredParam<SubdomainName>("primary_subdomain", "The name of the primary subdomain.");
53  params.addRequiredParam<SubdomainName>("secondary_subdomain",
54  "The name of the secondary subdomain.");
55  params.addParam<bool>(
56  "periodic",
57  false,
58  "Whether this constraint is going to be used to enforce a periodic condition. This has the "
59  "effect of changing the normals vector for projection from outward to inward facing");
60 
61  params.addParam<bool>(
62  "debug_mesh",
63  false,
64  "Whether this constraint is going to enable mortar segment mesh debug information. An exodus"
65  "file will be generated if the user sets this flag to true");
66 
67  params.addParam<bool>(
68  "correct_edge_dropping",
69  false,
70  "Whether to enable correct edge dropping treatment for mortar constraints. When disabled "
71  "any Lagrange Multiplier degree of freedom on a secondary element without full primary "
72  "contributions will be set (strongly) to 0.");
73 
74  params.addParam<bool>(
75  "interpolate_normals",
76  true,
77  "Whether to interpolate the nodal normals (e.g. classic idea of evaluating field at "
78  "quadrature points). If this is set to false, then non-interpolated nodal normals will be "
79  "used, and then the _normals member should be indexed with _i instead of _qp");
80 
81  params.addParam<bool>("ghost_point_neighbors",
82  false,
83  "Whether we should ghost point neighbors of secondary face elements, and "
84  "consequently also their mortar interface couples.");
85  params.addParam<Real>(
86  "minimum_projection_angle",
87  40.0,
88  "Parameter to control which angle (in degrees) is admissible for the creation of mortar "
89  "segments. If set to a value close to zero, very oblique projections are allowed, which "
90  "can result in mortar segments solving physics not meaningfully, and overprojection of "
91  "primary nodes onto the mortar segment mesh in extreme cases. This parameter is mostly "
92  "intended for mortar mesh debugging purposes in two dimensions.");
93 
94  params.addParam<bool>(
95  "ghost_higher_d_neighbors",
96  false,
97  "Whether we should ghost higher-dimensional neighbors. This is necessary when we are doing "
98  "second order mortar with finite volume primal variables, because in order for the method to "
99  "be second order we must use cell gradients, which couples in the neighbor cells.");
100 
101  return params;
102 }
103 
104 // Standard constructor
106  : _mci_fe_problem(*moose_object->getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
107  _mci_subproblem(*moose_object->getCheckedPointerParam<SubProblem *>("_subproblem")),
108  _mci_tid(moose_object->getParam<THREAD_ID>("_tid")),
109  _mci_mesh(_mci_subproblem.mesh()),
110  // all geometric assembly information should be correct for nl system number 0
111  _mci_assembly(_mci_subproblem.assembly(_mci_tid, 0)),
112  _mortar_data(_mci_fe_problem.mortarData()),
113  _secondary_id(
114  _mci_mesh.getBoundaryID(moose_object->getParam<BoundaryName>("secondary_boundary"))),
115  _primary_id(_mci_mesh.getBoundaryID(moose_object->getParam<BoundaryName>("primary_boundary"))),
116  _secondary_subdomain_id(
117  _mci_mesh.getSubdomainID(moose_object->getParam<SubdomainName>("secondary_subdomain"))),
118  _primary_subdomain_id(
119  _mci_mesh.getSubdomainID(moose_object->getParam<SubdomainName>("primary_subdomain"))),
120  _secondary_set({_secondary_id}),
121  _interpolate_normals(moose_object->getParam<bool>("interpolate_normals")),
122  _phys_points_secondary(_mci_assembly.qPointsFace()),
123  _phys_points_primary(_mci_assembly.qPointsFaceNeighbor()),
124  _qrule_msm(_mci_assembly.qRuleMortar()),
125  _qrule_face(_mci_assembly.qRuleFace()),
126  _lower_secondary_elem(_mci_assembly.lowerDElem()),
127  _lower_primary_elem(_mci_assembly.neighborLowerDElem()),
128  _JxW_msm(_mci_assembly.jxWMortar()),
129  _msm_elem(_mci_assembly.msmElem())
130 {
131  const bool displaced = moose_object->isParamValid("use_displaced_mesh")
132  ? moose_object->getParam<bool>("use_displaced_mesh")
133  : false;
134 
135  // Create the mortar interface if it hasn't already been created
136  _mci_fe_problem.createMortarInterface(
137  std::make_pair(_primary_id, _secondary_id),
138  std::make_pair(_primary_subdomain_id, _secondary_subdomain_id),
139  displaced,
140  moose_object->getParam<bool>("periodic"),
141  moose_object->getParam<bool>("debug_mesh"),
142  moose_object->getParam<bool>("correct_edge_dropping"),
143  moose_object->getParam<Real>("minimum_projection_angle"));
144 
145  _amg = &_mci_fe_problem.getMortarInterface(
146  std::make_pair(_primary_id, _secondary_id),
147  std::make_pair(_primary_subdomain_id, _secondary_subdomain_id),
148  displaced);
149 
150  const auto & secondary_set = _mortar_data.getHigherDimSubdomainIDs(_secondary_subdomain_id);
151  const auto & primary_set = _mortar_data.getHigherDimSubdomainIDs(_primary_subdomain_id);
152 
153  std::set_union(secondary_set.begin(),
154  secondary_set.end(),
155  primary_set.begin(),
156  primary_set.end(),
157  std::inserter(_higher_dim_subdomain_ids, _higher_dim_subdomain_ids.begin()));
158  _boundary_ids = {_secondary_id, _primary_id};
159 }
160 
161 void
163 {
164  if (interpolateNormals())
166  else
168 }
169 
170 void
172  ADReal & dual_number)
173 {
174  auto md_it = dual_number.derivatives().nude_data().begin();
175  auto mi_it = dual_number.derivatives().nude_indices().begin();
176 
177  auto d_it = dual_number.derivatives().nude_data().begin();
178 
179  for (auto i_it = dual_number.derivatives().nude_indices().begin();
180  i_it != dual_number.derivatives().nude_indices().end();
181  ++i_it, ++d_it)
182  if (*i_it != remove_derivative_index)
183  {
184  *mi_it = *i_it;
185  *md_it = *d_it;
186  ++mi_it;
187  ++md_it;
188  }
189 
190  std::size_t n_indices = md_it - dual_number.derivatives().nude_data().begin();
191  dual_number.derivatives().nude_indices().resize(n_indices);
192  dual_number.derivatives().nude_data().resize(n_indices);
193 }
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:47
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:209
uint8_t dof_id_type