Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
GradientJumpIndicator.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 "GradientJumpIndicator.h"
11 
13 
16 {
18  params.addClassDescription(
19  "Compute the jump of the solution gradient across element boundaries.");
20  params.addParam<bool>("variable_is_FV",
21  false,
22  "Whether the solution variable is using a finite volume discretization");
23 
24  // We need more ghosting to compute finite volume gradients across from a boundary
25  // We do not use skewness correction here, therefore avoiding needing three layers
27  "ElementSideNeighborLayers",
29  [](const InputParameters & obj_params, InputParameters & rm_params) {
30  rm_params.set<unsigned short>("layers") = obj_params.get<bool>("variable_is_FV") ? 2 : 1;
31  });
32 
33  return params;
34 }
35 
37  : InternalSideIndicator(parameters)
38 {
39 }
40 
41 Real
43 {
44  Real jump = 0;
45  // If the variable is not defined in the neighbor cell, we cant define the block
46  // If the indicator is not defined in the neighbor cell, we should not be looking at it
47  mooseAssert(_neighbor_elem, "Should have a neighbor");
48  if (_var.hasBlocks(_neighbor_elem->subdomain_id()) && hasBlocks(_neighbor_elem->subdomain_id()))
49  {
50  if (_var.isFV())
51  jump =
53  Moose::ElemArg{_current_elem, /*correct_skewness=*/false}, Moose::currentState())) -
55  _var.gradient(Moose::ElemArg{_neighbor_elem, false}, Moose::currentState()))) *
56  _normals[_qp];
57  else
58  jump = (_grad_u[_qp] - _grad_u_neighbor[_qp]) * _normals[_qp];
59  }
60 
61  return jump * jump;
62 }
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.
const OutputTools< ComputeValueType >::VariableGradient & _grad_u_neighbor
Holds the current solution gradient at the current quadrature point.
virtual bool isFV() const
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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.
static InputParameters validParams()
Factory constructor initializes all internal references needed for indicator computation.
The InternalSideIndicator class is responsible for calculating the residuals for various physics on i...
A structure that is used to evaluate Moose functors logically at an element/cell center.
GradientJumpIndicator(const InputParameters &parameters)
unsigned int _qp
The current quadrature point.
MooseVariableField< ComputeValueType > & _var
const Elem *const & _neighbor_elem
The neighbor element across from the current side.
const OutputTools< ComputeValueType >::VariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
registerMooseObject("MooseApp", GradientJumpIndicator)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool hasBlocks(const SubdomainID id) const override
Returns whether the functor is defined on this block.
GradientType gradient(const ElemArg &elem, const StateArg &state) const
Same as their evaluateGradient overloads with the same arguments but allows for caching implementatio...
Definition: MooseFunctor.h:847
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
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...
static InputParameters validParams()
bool hasBlocks(const SubdomainName &name) const
Test if the supplied block name is valid for this object.
virtual Real computeQpIntegral() override
The virtual function you will want to override to compute error contributions.
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
StateArg currentState()