https://mooseframework.inl.gov
Function.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 "Function.h"
11 #include "FEProblemBase.h"
12 
13 using namespace Moose;
14 
17 {
19  params += SetupInterface::validParams();
20 
21  // Functions should be executed on the fly
22  params.suppressParameter<ExecFlagEnum>("execute_on");
23  params.registerBase("Function");
24 
25  return params;
26 }
27 
29  : MooseObject(parameters),
30  SetupInterface(this),
31  TransientInterface(this),
33  UserObjectInterface(this),
34  Restartable(this, "Functions"),
35  MeshChangedInterface(parameters),
36  ScalarCoupleable(this),
38 {
39 }
40 
42 
43 Real
44 Function::value(Real /*t*/, const Point & /*p*/) const
45 {
46  mooseError("value method not implemented");
47  return 0.0;
48 }
49 
50 ADReal
51 Function::value(const ADReal & t, const ADPoint & p) const
52 {
53  const auto rt = MetaPhysicL::raw_value(t);
54  const auto rp = MetaPhysicL::raw_value(p);
55  const auto grad = gradient(rt, rp);
56  ADReal ret = value(rt, rp);
57  ret.derivatives() = grad(0) * p(0).derivatives()
58 #if LIBMESH_DIM > 1
59  + grad(1) * p(1).derivatives()
60 #endif
61 #if LIBMESH_DIM > 2
62  + grad(2) * p(2).derivatives()
63 #endif
64  + timeDerivative(rt, rp) * t.derivatives();
65  return ret;
66 }
67 
69 Function::value(const ChainedReal & t) const
70 {
71  static const Point p;
72  return {value(t.value(), p), timeDerivative(t.value(), p) * t.derivatives()};
73 }
74 
76 Function::gradient(Real /*t*/, const Point & /*p*/) const
77 {
78  mooseError("gradient method not implemented");
79  return RealGradient(0, 0, 0);
80 }
81 
82 Real
83 Function::timeDerivative(Real /*t*/, const Point & /*p*/) const
84 {
85  mooseError("timeDerivative method not implemented");
86  return 0;
87 }
88 
89 Real
90 Function::timeIntegral(Real /*t1*/, Real /*t2*/, const Point & /*p*/) const
91 {
92  mooseError("timeIntegral() not implemented.");
93 }
94 
96 Function::vectorValue(Real /*t*/, const Point & /*p*/) const
97 {
98  mooseError("vectorValue method not implemented");
99  return RealVectorValue(0, 0, 0);
100 }
101 
103 Function::curl(Real /*t*/, const Point & /*p*/) const
104 {
105  mooseError("curl method not implemented");
106  return RealVectorValue(0, 0, 0);
107 }
108 
109 Real
110 Function::div(Real /*t*/, const Point & /*p*/) const
111 {
112  mooseError("div method not implemented");
113  return 0;
114 }
115 
116 Real
118 {
119  mooseError("Integral method not implemented for function ", name());
120  return 0;
121 }
122 
123 Real
125 {
126  mooseError("Average method not implemented for function ", name());
127  return 0;
128 }
129 
130 template <typename R>
131 typename Function::ValueType
132 Function::evaluateHelper(const R & r, const Moose::StateArg & state) const
133 {
134  return value(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
135 }
136 
137 typename Function::ValueType
138 Function::evaluate(const ElemArg & elem_arg, const Moose::StateArg & state) const
139 {
140  return evaluateHelper(elem_arg, state);
141 }
142 
143 typename Function::ValueType
144 Function::evaluate(const FaceArg & face, const Moose::StateArg & state) const
145 {
146  if (face.face_side && face.fi->neighborPtr() &&
147  (face.fi->elem().subdomain_id() != face.fi->neighbor().subdomain_id()))
148  {
149  // Some users like to put discontinuities in their functions at subdomain changes in which case
150  // in order to always get the proper discontinuous effect we should evaluate ever so slightly
151  // off the face. Consider evaluation of: if(x < 0, -1, 1) if the face centroid is right at x ==
152  // 0 for example. The user likely doesn't want you to return 1 if they've asked for a 0-
153  // evaluation
154  //
155  // I can't quite tell but I think the tolerance for comparing coordinates (x, y, z, t) in
156  // fparser is ~1e-9 so we need to use something larger than that. The comparison is absolute
157  static constexpr Real offset_tolerance = 1e-8;
158  auto offset = offset_tolerance * face.fi->normal();
159  if (face.face_side == face.fi->elemPtr())
160  offset *= -1;
161  return value(_ti_feproblem.getTimeFromStateArg(state), face.getPoint() + offset);
162  }
163  else
164  return value(_ti_feproblem.getTimeFromStateArg(state), face.getPoint());
165 }
166 
167 typename Function::ValueType
168 Function::evaluate(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
169 {
170  return evaluateHelper(elem_qp, state);
171 }
172 
173 typename Function::ValueType
174 Function::evaluate(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
175 {
176  return evaluateHelper(elem_side_qp, state);
177 }
178 
179 typename Function::ValueType
180 Function::evaluate(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
181 {
182  return evaluateHelper(elem_point_arg, state);
183 }
184 
185 typename Function::ValueType
186 Function::evaluate(const NodeArg & node_arg, const Moose::StateArg & state) const
187 {
188  return evaluateHelper(node_arg, state);
189 }
190 
191 template <typename R>
192 typename Function::GradientType
193 Function::evaluateGradientHelper(const R & r, const Moose::StateArg & state) const
194 {
195  return gradient(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
196 }
197 
198 typename Function::GradientType
199 Function::evaluateGradient(const ElemArg & elem_arg, const Moose::StateArg & state) const
200 {
201  return evaluateGradientHelper(elem_arg, state);
202 }
203 
204 typename Function::GradientType
205 Function::evaluateGradient(const FaceArg & face, const Moose::StateArg & state) const
206 {
207  return evaluateGradientHelper(face, state);
208 }
209 
210 typename Function::GradientType
211 Function::evaluateGradient(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
212 {
213  return evaluateGradientHelper(elem_qp, state);
214 }
215 
216 typename Function::GradientType
217 Function::evaluateGradient(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
218 {
219  return evaluateGradientHelper(elem_side_qp, state);
220 }
221 
222 typename Function::GradientType
223 Function::evaluateGradient(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
224 {
225  return evaluateGradientHelper(elem_point_arg, state);
226 }
227 
228 typename Function::GradientType
229 Function::evaluateGradient(const NodeArg & node_arg, const Moose::StateArg & state) const
230 {
231  return evaluateGradientHelper(node_arg, state);
232 }
233 
234 template <typename R>
235 typename Function::DotType
236 Function::evaluateDotHelper(const R & r, const Moose::StateArg & state) const
237 {
238  return timeDerivative(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
239 }
240 
241 typename Function::DotType
242 Function::evaluateDot(const ElemArg & elem_arg, const Moose::StateArg & state) const
243 {
244  return evaluateDotHelper(elem_arg, state);
245 }
246 
247 typename Function::DotType
248 Function::evaluateDot(const FaceArg & face, const Moose::StateArg & state) const
249 {
250  return evaluateDotHelper(face, state);
251 }
252 
253 typename Function::DotType
254 Function::evaluateDot(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
255 {
256  return evaluateDotHelper(elem_qp, state);
257 }
258 
259 typename Function::DotType
260 Function::evaluateDot(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
261 {
262  return evaluateDotHelper(elem_side_qp, state);
263 }
264 
265 typename Function::DotType
266 Function::evaluateDot(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
267 {
268  return evaluateDotHelper(elem_point_arg, state);
269 }
270 
271 typename Function::DotType
272 Function::evaluateDot(const NodeArg & node_arg, const Moose::StateArg & state) const
273 {
274  return evaluateDotHelper(node_arg, state);
275 }
276 
277 void
279 {
281 }
282 
283 void
285 {
287 }
288 
289 void
291 {
293 }
294 
295 void
297 {
299 }
std::string name(const ElemQuality q)
A MultiMooseEnum object to hold "execute_on" flags.
Definition: ExecFlagEnum.h:21
A class for creating restricted objects.
Definition: Restartable.h:28
FEProblemBase & _ti_feproblem
Base class template for functor objects.
Definition: MooseFunctor.h:137
virtual Real timeIntegral(Real t1, Real t2, const Point &p) const
Computes the time integral at a spatial point between two time values.
Definition: Function.C:90
virtual Real div(Real t, const Point &p) const
Override this to evaluate the divergence of the vector function at a point (t,x,y,z), by default this returns zero, you must override it.
Definition: Function.C:110
const libMesh::Elem * face_side
A member that can be used to indicate whether there is a sidedness to this face.
DotType evaluateDotHelper(const R &r, const Moose::StateArg &state) const
Definition: Function.C:236
virtual RealVectorValue curl(Real t, const Point &p) const
Override this to evaluate the curl of the vector function at a point (t,x,y,z), by default this retur...
Definition: Function.C:103
virtual Real timeDerivative(Real t, const Point &p) const
Get the time derivative of the function.
Definition: Function.C:83
const Elem & elem() const
Definition: FaceInfo.h:81
DualNumber< Real, Real > ChainedReal
Definition: ChainedReal.h:30
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...
A structure that is used to evaluate Moose functors at an arbitrary physical point contained within a...
void residualSetup() override final
Gets called just before the residual is computed and before this object is asked to do its job...
Definition: Function.C:284
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
virtual const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:57
typename FunctorReturnType< Real, FunctorEvaluationKind::Gradient >::type GradientType
This rigmarole makes it so that a user can create functors that return containers (std::vector...
Definition: MooseFunctor.h:149
Function(const InputParameters &parameters)
Definition: Function.C:28
virtual Real average() const
Returns the average of the function over its domain.
Definition: Function.C:124
void suppressParameter(const std::string &name)
This method suppresses an inherited parameter so that it isn&#39;t required or valid in the derived class...
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
Interface for objects that needs transient capabilities.
const Elem * neighborPtr() const
Definition: FaceInfo.h:84
Interface for notifications that the mesh has changed.
A structure defining a "face" evaluation calling argument for Moose functors.
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
Real getTimeFromStateArg(const Moose::StateArg &state) const
Returns the time associated with the requested state.
const FaceInfo * fi
a face information object which defines our location in space
libMesh::Point getPoint() const
const Elem & neighbor() const
Definition: FaceInfo.h:216
A structure that is used to evaluate Moose functors logically at an element/cell center.
Argument for requesting functor evaluation at a quadrature point location in an element.
Interface for objects that need to use UserObjects.
const Point & normal() const
Returns the unit normal vector for the face oriented outward from the face&#39;s elem element...
Definition: FaceInfo.h:68
void jacobianSetup() override final
Gets called just before the Jacobian is computed and before this object is asked to do its job...
Definition: Function.C:290
void customSetup(const ExecFlagType &exec_type) override final
Gets called in FEProblemBase::execute() for execute flags other than initial, timestep_begin, nonlinear, linear and subdomain.
Definition: Function.C:296
void timestepSetup() override
Gets called at the beginning of the timestep before this object is asked to do its job...
Definition: Function.C:278
static InputParameters validParams()
const Elem * elemPtr() const
Definition: FaceInfo.h:82
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
GradientType evaluateGradientHelper(const R &r, const Moose::StateArg &state) const
Definition: Function.C:193
ValueType evaluate(const ElemArg &elem, const Moose::StateArg &state) const override final
Evaluate the functor with a given element.
Definition: Function.C:138
Class for containing MooseEnum item information.
Definition: MooseEnumItem.h:18
virtual RealGradient gradient(Real t, const Point &p) const
Function objects can optionally provide a gradient at a point.
Definition: Function.C:76
DotType evaluateDot(const ElemArg &elem, const Moose::StateArg &state) const override final
Evaluate the functor time derivative with a given element.
Definition: Function.C:242
Interface for objects that needs scalar coupling capabilities.
virtual RealVectorValue vectorValue(Real t, const Point &p) const
Override this to evaluate the vector function at a point (t,x,y,z), by default this returns a zero ve...
Definition: Function.C:96
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
virtual Real integral() const
Returns the integral of the function over its domain.
Definition: Function.C:117
State argument for evaluating functors.
virtual ~Function()
Function destructor.
Definition: Function.C:41
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
static InputParameters validParams()
Definition: MooseObject.C:25
virtual Real value(Real t, const Point &p) const
Override this to evaluate the scalar function at point (t,x,y,z), by default this returns zero...
Definition: Function.C:44
GradientType evaluateGradient(const ElemArg &elem, const Moose::StateArg &state) const override final
Evaluate the functor gradient with a given element.
Definition: Function.C:199
static InputParameters validParams()
Class constructor.
Definition: Function.C:16
Argument for requesting functor evaluation at quadrature point locations on an element side...
ValueType evaluateHelper(const R &r, const Moose::StateArg &state) const
Definition: Function.C:132
Interface class for classes which interact with Postprocessors.