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 : #pragma once 11 : 12 : #include "MooseTypes.h" 13 : 14 : #include <set> 15 : 16 : class InputParameters; 17 : class MooseObject; 18 : class MooseMesh; 19 : class MortarData; 20 : class AutomaticMortarGeneration; 21 : class SubProblem; 22 : class Assembly; 23 : 24 : namespace libMesh 25 : { 26 : class QBase; 27 : } 28 : 29 : /** 30 : * An interface for accessing mortar mesh data 31 : * 32 : * \note mci is shorthand for mortar consumer interface, so \p _mci_problem indicates 33 : * the mortar consumer interface's problem 34 : */ 35 : class MortarConsumerInterface 36 : { 37 : public: 38 : MortarConsumerInterface(const MooseObject * moose_object); 39 : 40 : static InputParameters validParams(); 41 : 42 : /** 43 : * @return The primary lower dimensional subdomain id 44 : */ 45 486 : SubdomainID primarySubdomain() const { return _primary_subdomain_id; } 46 : 47 : /** 48 : * @return The secondary lower dimensional subdomain id 49 : */ 50 486 : SubdomainID secondarySubdomain() const { return _secondary_subdomain_id; } 51 : 52 : /** 53 : * @return Whether this object lies on the given primary-secondary boundary pair 54 : */ 55 : bool onInterface(BoundaryID primary_boundary_id, BoundaryID secondary_boundary_id) const; 56 : 57 : protected: 58 : const std::set<SubdomainID> & getHigherDimSubdomainIDs() const 59 : { 60 : return _higher_dim_subdomain_ids; 61 : } 62 1232 : const std::set<BoundaryID> & getBoundaryIDs() const { return _boundary_ids; } 63 : 64 : /** 65 : * Retrieve the automatic mortar generation object associated with this constraint 66 : */ 67 : const AutomaticMortarGeneration & amg() const; 68 : 69 : /** 70 : * Whether to interpolate the nodal normals (e.g. classic idea of evaluating field at quadrature 71 : * points). If this is set to false, then non-interpolated nodal normals will be used, and then 72 : * the _normals member should be indexed with _i instead of _qp 73 : */ 74 279422 : bool interpolateNormals() const { return _interpolate_normals; } 75 : 76 : /** 77 : * Get rid of AD derivative entries by dof index 78 : */ 79 : static void trimDerivative(dof_id_type remove_derivative_index, ADReal & dual_number); 80 : 81 : /** 82 : * Get rid of interior node variable's derivatives 83 : */ 84 : template <typename Variables, typename DualNumbers> 85 : static void 86 : trimInteriorNodeDerivatives(const std::map<unsigned int, unsigned int> & primary_ip_lowerd_map, 87 : const Variables & moose_var, 88 : DualNumbers & ad_vars, 89 : const bool is_secondary); 90 : 91 : FEProblemBase & _mci_fe_problem; 92 : SubProblem & _mci_subproblem; 93 : const THREAD_ID _mci_tid; 94 : 95 : /// Mesh to query for boundary and subdomain ID information 96 : MooseMesh & _mci_mesh; 97 : 98 : Assembly & _mci_assembly; 99 : 100 : /// A reference to the mortar data object that holds all the mortar 101 : /// mesh information 102 : const MortarData & _mortar_data; 103 : 104 : /// Boundary ID for the secondary surface 105 : const BoundaryID _secondary_id; 106 : 107 : /// Boundary ID for the primary surface 108 : const BoundaryID _primary_id; 109 : 110 : /// Subdomain ID for the secondary surface 111 : const SubdomainID _secondary_subdomain_id; 112 : 113 : /// Subdomain ID for the primary surface 114 : const SubdomainID _primary_subdomain_id; 115 : 116 : /// the secondaryid set 117 : std::set<BoundaryID> _secondary_set; 118 : 119 : /// the union of the secondary and primary boundary ids 120 : std::set<BoundaryID> _boundary_ids; 121 : 122 : /// the higher dimensional subdomain ids corresponding to the interior parents 123 : std::set<SubdomainID> _higher_dim_subdomain_ids; 124 : 125 : /// Whether to interpolate the nodal normals 126 : const bool _interpolate_normals; 127 : 128 : /// The locations of the quadrature points on the interior secondary elements 129 : const MooseArray<Point> & _phys_points_secondary; 130 : 131 : /// The locations of the quadrature points on the interior primary elements 132 : const MooseArray<Point> & _phys_points_primary; 133 : 134 : /// The quadrature rule on the mortar segment element 135 : const libMesh::QBase * const & _qrule_msm; 136 : 137 : /// The arbitrary quadrature rule on the lower dimensional secondary face 138 : const libMesh::QBase * const & _qrule_face; 139 : 140 : /// The secondary face lower dimensional element (not the mortar element!). The mortar element 141 : /// lives on the secondary side of the mortar interface and *may* correspond to \p 142 : /// _lower_secondary_elem under the very specific circumstance that the nodes on the primary side 143 : /// of the mortar interface exactly project onto the secondary side of the mortar interface. In 144 : /// general projection of primary nodes will split the face elements on the secondary side of the 145 : /// interface. It is these split elements that are the mortar segment mesh elements 146 : Elem const * const & _lower_secondary_elem; 147 : 148 : /// The primary face lower dimensional element (not the mortar element!). The mortar element 149 : /// lives on the secondary side of the mortar interface and *may* correspond to \p 150 : /// _lower_secondary_elem under the very specific circumstance that the nodes on the primary side 151 : /// of the mortar interface exactly project onto the secondary side of the mortar interface. In 152 : /// general projection of primary nodes will split the face elements on the secondary side of the 153 : /// interface. It is these split elements that are the mortar segment mesh elements 154 : Elem const * const & _lower_primary_elem; 155 : 156 : /// The element Jacobian times weights 157 : const std::vector<Real> & _JxW_msm; 158 : 159 : /// The current mortar segment element 160 : const Elem * const & _msm_elem; 161 : 162 : /// the normals 163 : std::vector<Point> _normals; 164 : 165 : private: 166 : /** 167 : * Set the normals vector 168 : */ 169 : void setNormals(); 170 : 171 : // Pointer to automatic mortar generation object to give constraints access to mortar geometry 172 : const AutomaticMortarGeneration * _amg; 173 : 174 : friend class ComputeMortarFunctor; 175 : friend class FEProblemBase; 176 : 177 : template <typename> 178 : friend class MortarNodalAuxKernelTempl; 179 : friend class MortarUserObjectThread; 180 : 181 : /** 182 : * Calls the reinitialization of mortar user objects 183 : * @see FEProblemBase::reinitMortarUserObjects 184 : */ 185 : friend void reinitMortarUserObjects(BoundaryID, BoundaryID, bool); 186 : }; 187 : 188 : inline const AutomaticMortarGeneration & 189 343090 : MortarConsumerInterface::amg() const 190 : { 191 : mooseAssert(_amg, "this should have been set in the constructor"); 192 343090 : return *_amg; 193 : } 194 : 195 : template <typename Variables, typename DualNumbers> 196 : void 197 41400 : MortarConsumerInterface::trimInteriorNodeDerivatives( 198 : const std::map<unsigned int, unsigned int> & domain_ip_lowerd_map, 199 : const Variables & moose_vars, 200 : DualNumbers & dual_numbers, 201 : const bool is_secondary) 202 : { 203 : // 204 : // Remove interior node variable's derivatives from AD objects. 205 : // 206 : 207 : mooseAssert(moose_vars.size(), "Should have passed at least one variable"); 208 41400 : const auto num_indices = is_secondary ? moose_vars[0]->dofIndices().size() 209 20700 : : moose_vars[0]->dofIndicesNeighbor().size(); 210 : #ifdef DEBUG 211 : for (const auto i : make_range(std::size_t(1), moose_vars.size())) 212 : if (auto * moose_var = moose_vars[i]) 213 : mooseAssert(is_secondary ? moose_var->dofIndices().size() 214 : : moose_var->dofIndicesNeighbor().size() == num_indices, 215 : "These must be the same for all passed in variables"); 216 : #endif 217 : 218 197160 : for (const auto dof_index : make_range(num_indices)) 219 155760 : if (!domain_ip_lowerd_map.count(dof_index)) 220 : { 221 218880 : for (const auto * const moose_var : moose_vars) 222 : { 223 : // It's valid for a user to pass a container that represents a LIBMESH_DIM vector of 224 : // component variables for which one or two of the variables may be null depending on the 225 : // mesh dimension in the simulation 226 145920 : if (!moose_var) 227 72960 : continue; 228 : 229 : mooseAssert(moose_var->isNodal(), 230 : "Trimming of interior node's derivatives is only supported for Lagrange " 231 : "elements in mortar objects"); 232 : 233 72960 : const auto remove_derivative_index = is_secondary 234 114360 : ? moose_var->dofIndices()[dof_index] 235 41400 : : moose_var->dofIndicesNeighbor()[dof_index]; 236 283200 : for (auto & dual_number : dual_numbers) 237 210240 : trimDerivative(remove_derivative_index, dual_number); 238 : } 239 : } 240 41400 : } 241 : 242 : inline bool 243 22 : MortarConsumerInterface::onInterface(const BoundaryID primary_boundary_id, 244 : const BoundaryID secondary_boundary_id) const 245 : { 246 22 : return (primary_boundary_id == _primary_id) && (secondary_boundary_id == _secondary_id); 247 : }