Line data Source code
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
11 : #include "MortarConsumerInterface.h"
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 :
22 : InputParameters
23 259700 : MortarConsumerInterface::validParams()
24 : {
25 : // Create InputParameters object that will be appended to the parameters for the inheriting object
26 259700 : InputParameters params = emptyInputParameters();
27 : // On a displaced mesh this will geometrically and algebraically ghost the entire interface
28 779100 : params.addRelationshipManager(
29 : "AugmentSparsityOnInterface",
30 : Moose::RelationshipManagerType::GEOMETRIC | Moose::RelationshipManagerType::ALGEBRAIC,
31 259700 : [](const InputParameters & obj_params, InputParameters & rm_params)
32 : {
33 8792 : rm_params.set<bool>("use_displaced_mesh") = obj_params.get<bool>("use_displaced_mesh");
34 8792 : rm_params.set<BoundaryName>("secondary_boundary") =
35 8792 : obj_params.get<BoundaryName>("secondary_boundary");
36 8792 : rm_params.set<BoundaryName>("primary_boundary") =
37 8792 : obj_params.get<BoundaryName>("primary_boundary");
38 8792 : rm_params.set<SubdomainName>("secondary_subdomain") =
39 8792 : obj_params.get<SubdomainName>("secondary_subdomain");
40 8792 : rm_params.set<SubdomainName>("primary_subdomain") =
41 8792 : obj_params.get<SubdomainName>("primary_subdomain");
42 4396 : rm_params.set<bool>("ghost_point_neighbors") =
43 4396 : obj_params.get<bool>("ghost_point_neighbors");
44 4396 : rm_params.set<bool>("ghost_higher_d_neighbors") =
45 4396 : obj_params.get<bool>("ghost_higher_d_neighbors");
46 4396 : });
47 :
48 1038800 : params.addRequiredParam<BoundaryName>("primary_boundary",
49 : "The name of the primary boundary sideset.");
50 1038800 : params.addRequiredParam<BoundaryName>("secondary_boundary",
51 : "The name of the secondary boundary sideset.");
52 1038800 : params.addRequiredParam<SubdomainName>("primary_subdomain", "The name of the primary subdomain.");
53 1038800 : params.addRequiredParam<SubdomainName>("secondary_subdomain",
54 : "The name of the secondary subdomain.");
55 779100 : params.addParam<bool>(
56 : "periodic",
57 519400 : 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 779100 : params.addParam<bool>(
62 : "debug_mesh",
63 519400 : 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 779100 : params.addParam<bool>(
68 : "correct_edge_dropping",
69 519400 : 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 779100 : params.addParam<bool>(
75 : "interpolate_normals",
76 519400 : 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 779100 : params.addParam<bool>("ghost_point_neighbors",
82 519400 : false,
83 : "Whether we should ghost point neighbors of secondary face elements, and "
84 : "consequently also their mortar interface couples.");
85 779100 : params.addParam<Real>(
86 : "minimum_projection_angle",
87 519400 : 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 519400 : params.addParam<bool>(
95 : "ghost_higher_d_neighbors",
96 519400 : 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 259700 : return params;
102 0 : }
103 :
104 : // Standard constructor
105 1448 : MortarConsumerInterface::MortarConsumerInterface(const MooseObject * moose_object)
106 5792 : : _mci_fe_problem(*moose_object->getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
107 5792 : _mci_subproblem(*moose_object->getCheckedPointerParam<SubProblem *>("_subproblem")),
108 2896 : _mci_tid(moose_object->getParam<THREAD_ID>("_tid")),
109 1448 : _mci_mesh(_mci_subproblem.mesh()),
110 : // all geometric assembly information should be correct for nl system number 0
111 1448 : _mci_assembly(_mci_subproblem.assembly(_mci_tid, 0)),
112 1448 : _mortar_data(_mci_fe_problem.mortarData()),
113 1448 : _secondary_id(
114 2896 : _mci_mesh.getBoundaryID(moose_object->getParam<BoundaryName>("secondary_boundary"))),
115 2896 : _primary_id(_mci_mesh.getBoundaryID(moose_object->getParam<BoundaryName>("primary_boundary"))),
116 1448 : _secondary_subdomain_id(
117 2896 : _mci_mesh.getSubdomainID(moose_object->getParam<SubdomainName>("secondary_subdomain"))),
118 1448 : _primary_subdomain_id(
119 2896 : _mci_mesh.getSubdomainID(moose_object->getParam<SubdomainName>("primary_subdomain"))),
120 2896 : _secondary_set({_secondary_id}),
121 2896 : _interpolate_normals(moose_object->getParam<bool>("interpolate_normals")),
122 1448 : _phys_points_secondary(_mci_assembly.qPointsFace()),
123 1448 : _phys_points_primary(_mci_assembly.qPointsFaceNeighbor()),
124 1448 : _qrule_msm(_mci_assembly.qRuleMortar()),
125 1448 : _qrule_face(_mci_assembly.qRuleFace()),
126 1448 : _lower_secondary_elem(_mci_assembly.lowerDElem()),
127 1448 : _lower_primary_elem(_mci_assembly.neighborLowerDElem()),
128 1448 : _JxW_msm(_mci_assembly.jxWMortar()),
129 2896 : _msm_elem(_mci_assembly.msmElem())
130 : {
131 2896 : const bool displaced = moose_object->isParamValid("use_displaced_mesh")
132 5792 : ? moose_object->getParam<bool>("use_displaced_mesh")
133 1448 : : false;
134 :
135 : // Create the mortar interface if it hasn't already been created
136 1448 : _mci_fe_problem.createMortarInterface(
137 1448 : std::make_pair(_primary_id, _secondary_id),
138 0 : std::make_pair(_primary_subdomain_id, _secondary_subdomain_id),
139 : displaced,
140 2896 : moose_object->getParam<bool>("periodic"),
141 2896 : moose_object->getParam<bool>("debug_mesh"),
142 2896 : moose_object->getParam<bool>("correct_edge_dropping"),
143 4344 : moose_object->getParam<Real>("minimum_projection_angle"));
144 :
145 1448 : _amg = &_mci_fe_problem.getMortarInterface(
146 1448 : std::make_pair(_primary_id, _secondary_id),
147 1448 : std::make_pair(_primary_subdomain_id, _secondary_subdomain_id),
148 : displaced);
149 :
150 1448 : const auto & secondary_set = _mortar_data.getHigherDimSubdomainIDs(_secondary_subdomain_id);
151 1448 : const auto & primary_set = _mortar_data.getHigherDimSubdomainIDs(_primary_subdomain_id);
152 :
153 2896 : std::set_union(secondary_set.begin(),
154 : secondary_set.end(),
155 : primary_set.begin(),
156 : primary_set.end(),
157 1448 : std::inserter(_higher_dim_subdomain_ids, _higher_dim_subdomain_ids.begin()));
158 1448 : _boundary_ids = {_secondary_id, _primary_id};
159 1448 : }
160 :
161 : void
162 301433 : MortarConsumerInterface::setNormals()
163 : {
164 301433 : if (interpolateNormals())
165 300077 : _normals = amg().getNormals(*_lower_secondary_elem, _qrule_face->get_points());
166 : else
167 1356 : _normals = amg().getNodalNormals(*_lower_secondary_elem);
168 301433 : }
169 :
170 : void
171 224728 : MortarConsumerInterface::trimDerivative(const dof_id_type remove_derivative_index,
172 : ADReal & dual_number)
173 : {
174 224728 : auto md_it = dual_number.derivatives().nude_data().begin();
175 224728 : auto mi_it = dual_number.derivatives().nude_indices().begin();
176 :
177 224728 : auto d_it = dual_number.derivatives().nude_data().begin();
178 :
179 224728 : for (auto i_it = dual_number.derivatives().nude_indices().begin();
180 633728 : i_it != dual_number.derivatives().nude_indices().end();
181 409000 : ++i_it, ++d_it)
182 409000 : if (*i_it != remove_derivative_index)
183 : {
184 381588 : *mi_it = *i_it;
185 381588 : *md_it = *d_it;
186 381588 : ++mi_it;
187 381588 : ++md_it;
188 : }
189 :
190 224728 : std::size_t n_indices = md_it - dual_number.derivatives().nude_data().begin();
191 224728 : dual_number.derivatives().nude_indices().resize(n_indices);
192 224728 : dual_number.derivatives().nude_data().resize(n_indices);
193 224728 : }
|