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 : #include "Function.h"
11 : #include "FEProblemBase.h"
12 :
13 : using namespace Moose;
14 :
15 : InputParameters
16 571358 : Function::validParams()
17 : {
18 571358 : InputParameters params = MooseObject::validParams();
19 571358 : params += SetupInterface::validParams();
20 :
21 : // Functions should be executed on the fly
22 571358 : params.suppressParameter<ExecFlagEnum>("execute_on");
23 571358 : params.registerBase("Function");
24 :
25 571358 : return params;
26 0 : }
27 :
28 52290 : Function::Function(const InputParameters & parameters)
29 : : MooseObject(parameters),
30 : SetupInterface(this),
31 : TransientInterface(this),
32 : PostprocessorInterface(this),
33 : UserObjectInterface(this),
34 : Restartable(this, "Functions"),
35 : MeshChangedInterface(parameters),
36 : ScalarCoupleable(this),
37 156870 : Moose::FunctorBase<Real>(name())
38 : {
39 104580 : }
40 :
41 49238 : Function::~Function() {}
42 :
43 : Real
44 0 : Function::value(Real /*t*/, const Point & /*p*/) const
45 : {
46 0 : mooseError("value method not implemented");
47 : return 0.0;
48 : }
49 :
50 : ADReal
51 4552800 : Function::value(const ADReal & t, const ADPoint & p) const
52 : {
53 4552800 : const auto rt = MetaPhysicL::raw_value(t);
54 4552800 : const auto rp = MetaPhysicL::raw_value(p);
55 4552800 : const auto grad = gradient(rt, rp);
56 4552800 : ADReal ret = value(rt, rp);
57 13658400 : ret.derivatives() = grad(0) * p(0).derivatives()
58 : #if LIBMESH_DIM > 1
59 18211200 : + grad(1) * p(1).derivatives()
60 : #endif
61 : #if LIBMESH_DIM > 2
62 18211200 : + grad(2) * p(2).derivatives()
63 : #endif
64 18211200 : + timeDerivative(rt, rp) * t.derivatives();
65 9105600 : return ret;
66 0 : }
67 :
68 : ChainedReal
69 0 : Function::value(const ChainedReal & t) const
70 : {
71 0 : static const Point p;
72 0 : return {value(t.value(), p), timeDerivative(t.value(), p) * t.derivatives()};
73 : }
74 :
75 : RealGradient
76 0 : Function::gradient(Real /*t*/, const Point & /*p*/) const
77 : {
78 0 : mooseError("gradient method not implemented");
79 : return RealGradient(0, 0, 0);
80 : }
81 :
82 : Real
83 0 : Function::timeDerivative(Real /*t*/, const Point & /*p*/) const
84 : {
85 0 : mooseError("timeDerivative method not implemented");
86 : return 0;
87 : }
88 :
89 : Real
90 0 : Function::timeIntegral(Real /*t1*/, Real /*t2*/, const Point & /*p*/) const
91 : {
92 0 : mooseError("timeIntegral() not implemented.");
93 : }
94 :
95 : RealVectorValue
96 0 : Function::vectorValue(Real /*t*/, const Point & /*p*/) const
97 : {
98 0 : mooseError("vectorValue method not implemented");
99 : return RealVectorValue(0, 0, 0);
100 : }
101 :
102 : RealVectorValue
103 0 : Function::curl(Real /*t*/, const Point & /*p*/) const
104 : {
105 0 : mooseError("curl method not implemented");
106 : return RealVectorValue(0, 0, 0);
107 : }
108 :
109 : Real
110 0 : Function::div(Real /*t*/, const Point & /*p*/) const
111 : {
112 0 : mooseError("div method not implemented");
113 : return 0;
114 : }
115 :
116 : Real
117 0 : Function::integral() const
118 : {
119 0 : mooseError("Integral method not implemented for function ", name());
120 : return 0;
121 : }
122 :
123 : Real
124 0 : Function::average() const
125 : {
126 0 : mooseError("Average method not implemented for function ", name());
127 : return 0;
128 : }
129 :
130 : template <typename R>
131 : typename Function::ValueType
132 5958225 : Function::evaluateHelper(const R & r, const Moose::StateArg & state) const
133 : {
134 5958225 : return value(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
135 : }
136 :
137 : typename Function::ValueType
138 2981587 : Function::evaluate(const ElemArg & elem_arg, const Moose::StateArg & state) const
139 : {
140 2981587 : return evaluateHelper(elem_arg, state);
141 : }
142 :
143 : typename Function::ValueType
144 1648071 : Function::evaluate(const FaceArg & face, const Moose::StateArg & state) const
145 : {
146 1671867 : if (face.face_side && face.fi->neighborPtr() &&
147 23796 : (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 23796 : auto offset = offset_tolerance * face.fi->normal();
159 23796 : if (face.face_side == face.fi->elemPtr())
160 11898 : offset *= -1;
161 23796 : return value(_ti_feproblem.getTimeFromStateArg(state), face.getPoint() + offset);
162 : }
163 : else
164 1624275 : return value(_ti_feproblem.getTimeFromStateArg(state), face.getPoint());
165 : }
166 :
167 : typename Function::ValueType
168 2291177 : Function::evaluate(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
169 : {
170 2291177 : return evaluateHelper(elem_qp, state);
171 : }
172 :
173 : typename Function::ValueType
174 533441 : Function::evaluate(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
175 : {
176 533441 : return evaluateHelper(elem_side_qp, state);
177 : }
178 :
179 : typename Function::ValueType
180 699 : Function::evaluate(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
181 : {
182 699 : return evaluateHelper(elem_point_arg, state);
183 : }
184 :
185 : typename Function::ValueType
186 151321 : Function::evaluate(const NodeArg & node_arg, const Moose::StateArg & state) const
187 : {
188 151321 : return evaluateHelper(node_arg, state);
189 : }
190 :
191 : template <typename R>
192 : typename Function::GradientType
193 6405 : Function::evaluateGradientHelper(const R & r, const Moose::StateArg & state) const
194 : {
195 6405 : return gradient(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
196 : }
197 :
198 : typename Function::GradientType
199 1 : Function::evaluateGradient(const ElemArg & elem_arg, const Moose::StateArg & state) const
200 : {
201 1 : return evaluateGradientHelper(elem_arg, state);
202 : }
203 :
204 : typename Function::GradientType
205 1 : Function::evaluateGradient(const FaceArg & face, const Moose::StateArg & state) const
206 : {
207 1 : return evaluateGradientHelper(face, state);
208 : }
209 :
210 : typename Function::GradientType
211 6401 : Function::evaluateGradient(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
212 : {
213 6401 : return evaluateGradientHelper(elem_qp, state);
214 : }
215 :
216 : typename Function::GradientType
217 1 : Function::evaluateGradient(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
218 : {
219 1 : return evaluateGradientHelper(elem_side_qp, state);
220 : }
221 :
222 : typename Function::GradientType
223 1 : Function::evaluateGradient(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
224 : {
225 1 : return evaluateGradientHelper(elem_point_arg, state);
226 : }
227 :
228 : typename Function::GradientType
229 0 : Function::evaluateGradient(const NodeArg & node_arg, const Moose::StateArg & state) const
230 : {
231 0 : return evaluateGradientHelper(node_arg, state);
232 : }
233 :
234 : template <typename R>
235 : typename Function::DotType
236 2565 : Function::evaluateDotHelper(const R & r, const Moose::StateArg & state) const
237 : {
238 2565 : return timeDerivative(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
239 : }
240 :
241 : typename Function::DotType
242 257 : Function::evaluateDot(const ElemArg & elem_arg, const Moose::StateArg & state) const
243 : {
244 257 : return evaluateDotHelper(elem_arg, state);
245 : }
246 :
247 : typename Function::DotType
248 1 : Function::evaluateDot(const FaceArg & face, const Moose::StateArg & state) const
249 : {
250 1 : return evaluateDotHelper(face, state);
251 : }
252 :
253 : typename Function::DotType
254 2305 : Function::evaluateDot(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
255 : {
256 2305 : return evaluateDotHelper(elem_qp, state);
257 : }
258 :
259 : typename Function::DotType
260 1 : Function::evaluateDot(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
261 : {
262 1 : return evaluateDotHelper(elem_side_qp, state);
263 : }
264 :
265 : typename Function::DotType
266 1 : Function::evaluateDot(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
267 : {
268 1 : return evaluateDotHelper(elem_point_arg, state);
269 : }
270 :
271 : typename Function::DotType
272 0 : Function::evaluateDot(const NodeArg & node_arg, const Moose::StateArg & state) const
273 : {
274 0 : return evaluateDotHelper(node_arg, state);
275 : }
276 :
277 : void
278 221751 : Function::timestepSetup()
279 : {
280 221751 : FunctorBase<Real>::timestepSetup();
281 221751 : }
282 :
283 : void
284 2789449 : Function::residualSetup()
285 : {
286 2789449 : FunctorBase<Real>::residualSetup();
287 2789449 : }
288 :
289 : void
290 413950 : Function::jacobianSetup()
291 : {
292 413950 : FunctorBase<Real>::jacobianSetup();
293 413950 : }
294 :
295 : void
296 1413070 : Function::customSetup(const ExecFlagType & exec_type)
297 : {
298 1413070 : FunctorBase<Real>::customSetup(exec_type);
299 1413070 : }
|