https://mooseframework.inl.gov
DGKernelBase.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 #include "DGKernelBase.h"
11 #include "MooseVariable.h"
12 #include "Problem.h"
13 #include "SubProblem.h"
14 #include "SystemBase.h"
15 #include "MaterialData.h"
16 #include "ParallelUniqueId.h"
17 
18 #include "libmesh/dof_map.h"
19 #include "libmesh/dense_vector.h"
20 #include "libmesh/numeric_vector.h"
21 #include "libmesh/dense_subvector.h"
22 #include "libmesh/libmesh_common.h"
23 #include "libmesh/quadrature.h"
24 
27 {
32  params.addParam<bool>("use_displaced_mesh",
33  false,
34  "Whether or not this object should use the "
35  "displaced mesh for computation. Note that in "
36  "the case this is true but no displacements "
37  "are provided in the Mesh block the "
38  "undisplaced mesh will still be used.");
39  params.addParamNamesToGroup("use_displaced_mesh", "Advanced");
40 
41  params.addParam<std::vector<AuxVariableName>>(
42  "save_in",
43  {},
44  "The name of auxiliary variables to save this Kernel's residual contributions to. "
45  " Everything about that variable must match everything about this variable (the "
46  "type, what blocks it's on, etc.)");
47  params.addParam<std::vector<AuxVariableName>>(
48  "diag_save_in",
49  {},
50  "The name of auxiliary variables to save this Kernel's diagonal Jacobian "
51  "contributions to. Everything about that variable must match everything "
52  "about this variable (the type, what blocks it's on, etc.)");
53  params.addParamNamesToGroup("diag_save_in save_in", "Advanced");
54 
55  // DG Kernels always need one layer of ghosting.
56  params.addRelationshipManager("ElementSideNeighborLayers",
60 
61  params.addParam<std::vector<BoundaryName>>(
62  "exclude_boundary", "The internal side sets to be excluded from this kernel.");
63  params.registerBase("DGKernel");
64  params.registerSystemAttributeName("DGKernel");
65 
66  return params;
67 }
68 
70  : NeighborResidualObject(parameters),
71  BlockRestrictable(this),
72  BoundaryRestrictable(this, false), // false for _not_ nodal
74  TwoMaterialPropertyInterface(this, blockIDs(), boundaryIDs()),
75  ElementIDInterface(this),
76  _current_elem(_assembly.elem()),
77  _current_elem_volume(_assembly.elemVolume()),
78 
79  _neighbor_elem(_assembly.neighbor()),
80  _neighbor_elem_volume(_assembly.neighborVolume()),
81 
82  _current_side(_assembly.side()),
83  _current_side_elem(_assembly.sideElem()),
84  _current_side_volume(_assembly.sideElemVolume()),
85 
86  _coord_sys(_assembly.coordSystem()),
87  _q_point(_assembly.qPointsFace()),
88  _qrule(_assembly.qRuleFace()),
89  _JxW(_assembly.JxWFace()),
90  _coord(_assembly.coordTransformation()),
91  _normals(_assembly.normals()),
92 
93  _save_in_strings(parameters.get<std::vector<AuxVariableName>>("save_in")),
94  _diag_save_in_strings(parameters.get<std::vector<AuxVariableName>>("diag_save_in"))
95 {
96  // Gather information on broken boundaries
97  std::vector<BoundaryName> bnd = isParamValid("exclude_boundary")
98  ? getParam<std::vector<BoundaryName>>("exclude_boundary")
99  : std::vector<BoundaryName>(0);
100  auto bnd_ids = _mesh.getBoundaryIDs(bnd);
101 
102  // check if the broken boundary ids are valid
103  auto & valid_ids = _mesh.meshSidesetIds();
104  std::vector<BoundaryName> diff;
105  for (unsigned int i = 0; i < bnd_ids.size(); ++i)
106  if (valid_ids.find(bnd_ids[i]) == valid_ids.end())
107  diff.push_back(bnd[i]);
108  if (!diff.empty())
109  {
110  auto msg = "DGKernel '" + name() +
111  "' contains the following boundary names that do not exist on the mesh: " +
112  Moose::stringify(diff, ", ");
113  paramError("exclude_boundary", msg);
114  }
115 
116  _excluded_boundaries.insert(bnd_ids.begin(), bnd_ids.end());
117 }
118 
119 void
121 {
122  if (!excludeBoundary())
123  {
125 
126  // Compute the residual for this element
128 
129  // Compute the residual for the neighbor
131  }
132 }
133 
134 void
136 {
137  if (!excludeBoundary())
138  {
140 
141  // Compute element-element Jacobian
143 
144  // Compute element-neighbor Jacobian
146 
147  // Compute neighbor-element Jacobian
149 
150  // Compute neighbor-neighbor Jacobian
152  }
153 }
154 
155 void
156 DGKernelBase::computeOffDiagJacobian(const unsigned int jvar_num)
157 {
158  if (!excludeBoundary())
159  {
160  const auto & jvar = getVariable(jvar_num);
161 
162  if (jvar_num == variable().number())
163  computeJacobian();
164  else
165  {
166  precalculateOffDiagJacobian(jvar_num);
167 
168  // Compute element-element Jacobian
170 
171  // Compute element-neighbor Jacobian
173 
174  // Compute neighbor-element Jacobian
176 
177  // Compute neighbor-neighbor Jacobian
179  }
180  }
181 }
182 
183 bool
185 {
186  if (_excluded_boundaries.empty())
187  return false;
188 
189  auto boundary_ids = _mesh.getBoundaryIDs(_current_elem, _current_side);
190  for (auto bid : boundary_ids)
191  if (_excluded_boundaries.find(bid) != _excluded_boundaries.end())
192  return true;
193 
194  // make sure we will also break on the neighboring side
195  unsigned int neighbor_side = _neighbor_elem->which_neighbor_am_i(_current_elem);
196  boundary_ids = _mesh.getBoundaryIDs(_neighbor_elem, neighbor_side);
197  for (auto bid : boundary_ids)
198  if (_excluded_boundaries.find(bid) != _excluded_boundaries.end())
199  return true;
200 
201  return false;
202 }
203 
204 void
205 DGKernelBase::prepareShapes(const unsigned int var_num)
206 {
208 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
Definition: DGKernelBase.C:156
std::set< BoundaryID > _excluded_boundaries
Broken boundaries.
Definition: DGKernelBase.h:145
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernelBase.C:26
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:439
static InputParameters validParams()
const std::set< BoundaryID > & meshSidesetIds() const
Returns a read-only reference to the set of sidesets currently present in the Mesh.
Definition: MooseMesh.C:3263
virtual void precalculateOffDiagJacobian(unsigned int)
virtual void precalculateResidual()
DGKernelBase(const InputParameters &parameters)
Definition: DGKernelBase.C:69
void registerSystemAttributeName(const std::string &value)
This method is used to define the MOOSE system name that is used by the TheWarehouse object for stori...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
/class BoundaryRestrictable /brief Provides functionality for limiting the object to certain boundary...
virtual void computeJacobian() override
Computes the jacobian for the current side.
Definition: DGKernelBase.C:135
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.
THREAD_ID _tid
The thread ID for this kernel.
virtual void computeElemNeighResidual(Moose::DGResidualType type)=0
Computes the residual for this element or the neighbor.
const MooseVariableFieldBase & getVariable(unsigned int jvar_num) const
Retrieve the variable object from our system associated with jvar_num.
static InputParameters validParams()
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:103
virtual void precalculateJacobian()
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
Definition: DGKernelBase.C:205
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
bool excludeBoundary() const
Check current element if it contains broken boundary.
Definition: DGKernelBase.C:184
std::string stringify(const T &t)
conversion to string
Definition: Conversion.h:64
static InputParameters validParams()
This is a base class for objects that can provide residual contributions for both local and neighbor ...
virtual void computeResidual() override
Computes the residual for the current side.
Definition: DGKernelBase.C:120
const Elem *const & _current_elem
Current element.
Definition: DGKernelBase.h:84
An interface that restricts an object to subdomains via the &#39;blocks&#39; input parameter.
const unsigned int & _current_side
Current side.
Definition: DGKernelBase.h:96
This interface is designed for DGKernel, InternalSideUserObject, InterfaceUserObject, where material properties on a side of both its primary side (face) and its secondary side (neighbor) all required.
const Elem *const & _neighbor_elem
The neighboring element.
Definition: DGKernelBase.h:90
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...
virtual const MooseVariableBase & variable() const =0
Returns the variable that this object operates on.
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:199
virtual void computeElemNeighJacobian(Moose::DGJacobianType type)=0
Computes the element/neighbor-element/neighbor Jacobian.
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const
Returns a vector of boundary IDs for the requested element on the requested side. ...
virtual void prepareFaceShapes(unsigned int var, const THREAD_ID tid)=0
const Elem & get(const ElemType type_in)
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...
virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar)=0
Computes the element-element off-diagonal Jacobian.