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 590147 : Function::validParams()
17 : {
18 590147 : InputParameters params = FunctionBase::validParams();
19 590147 : return params;
20 : }
21 :
22 57400 : Function::Function(const InputParameters & parameters)
23 229600 : : FunctionBase(parameters), Moose::FunctorBase<Real>(name())
24 : {
25 114800 : }
26 :
27 54085 : Function::~Function() {}
28 :
29 : Real
30 0 : Function::value(Real /*t*/, const Point & /*p*/) const
31 : {
32 0 : mooseError("value method not implemented");
33 : return 0.0;
34 : }
35 :
36 : ADReal
37 4552800 : Function::value(const ADReal & t, const ADPoint & p) const
38 : {
39 4552800 : const auto rt = MetaPhysicL::raw_value(t);
40 4552800 : const auto rp = MetaPhysicL::raw_value(p);
41 4552800 : const auto grad = gradient(rt, rp);
42 4552800 : ADReal ret = value(rt, rp);
43 13658400 : ret.derivatives() = grad(0) * p(0).derivatives()
44 : #if LIBMESH_DIM > 1
45 18211200 : + grad(1) * p(1).derivatives()
46 : #endif
47 : #if LIBMESH_DIM > 2
48 18211200 : + grad(2) * p(2).derivatives()
49 : #endif
50 18211200 : + timeDerivative(rt, rp) * t.derivatives();
51 9105600 : return ret;
52 0 : }
53 :
54 : ChainedReal
55 0 : Function::value(const ChainedReal & t) const
56 : {
57 0 : static const Point p;
58 0 : return {value(t.value(), p), timeDerivative(t.value(), p) * t.derivatives()};
59 : }
60 :
61 : RealGradient
62 0 : Function::gradient(Real /*t*/, const Point & /*p*/) const
63 : {
64 0 : mooseError("gradient method not implemented");
65 : return RealGradient(0, 0, 0);
66 : }
67 :
68 : Real
69 0 : Function::timeDerivative(Real /*t*/, const Point & /*p*/) const
70 : {
71 0 : mooseError("timeDerivative method not implemented");
72 : return 0;
73 : }
74 :
75 : Real
76 0 : Function::timeIntegral(Real /*t1*/, Real /*t2*/, const Point & /*p*/) const
77 : {
78 0 : mooseError("timeIntegral() not implemented.");
79 : }
80 :
81 : RealVectorValue
82 0 : Function::vectorValue(Real /*t*/, const Point & /*p*/) const
83 : {
84 0 : mooseError("vectorValue method not implemented");
85 : return RealVectorValue(0, 0, 0);
86 : }
87 :
88 : RealVectorValue
89 0 : Function::curl(Real /*t*/, const Point & /*p*/) const
90 : {
91 0 : mooseError("curl method not implemented");
92 : return RealVectorValue(0, 0, 0);
93 : }
94 :
95 : Real
96 0 : Function::div(Real /*t*/, const Point & /*p*/) const
97 : {
98 0 : mooseError("div method not implemented");
99 : return 0;
100 : }
101 :
102 : Real
103 0 : Function::integral() const
104 : {
105 0 : mooseError("Integral method not implemented for function ", name());
106 : return 0;
107 : }
108 :
109 : Real
110 0 : Function::average() const
111 : {
112 0 : mooseError("Average method not implemented for function ", name());
113 : return 0;
114 : }
115 :
116 : template <typename R>
117 : typename Function::ValueType
118 6242419 : Function::evaluateHelper(const R & r, const Moose::StateArg & state) const
119 : {
120 6242419 : return value(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
121 : }
122 :
123 : typename Function::ValueType
124 3015641 : Function::evaluate(const ElemArg & elem_arg, const Moose::StateArg & state) const
125 : {
126 3015641 : return evaluateHelper(elem_arg, state);
127 : }
128 :
129 : typename Function::ValueType
130 1711444 : Function::evaluate(const FaceArg & face, const Moose::StateArg & state) const
131 : {
132 1735936 : if (face.face_side && face.fi->neighborPtr() &&
133 24492 : (face.fi->elem().subdomain_id() != face.fi->neighbor().subdomain_id()))
134 : {
135 : // Some users like to put discontinuities in their functions at subdomain changes in which case
136 : // in order to always get the proper discontinuous effect we should evaluate ever so slightly
137 : // off the face. Consider evaluation of: if(x < 0, -1, 1) if the face centroid is right at x ==
138 : // 0 for example. The user likely doesn't want you to return 1 if they've asked for a 0-
139 : // evaluation
140 : //
141 : // I can't quite tell but I think the tolerance for comparing coordinates (x, y, z, t) in
142 : // fparser is ~1e-9 so we need to use something larger than that. The comparison is absolute
143 : static constexpr Real offset_tolerance = 1e-8;
144 24492 : auto offset = offset_tolerance * face.fi->normal();
145 24492 : if (face.face_side == face.fi->elemPtr())
146 12246 : offset *= -1;
147 24492 : return value(_ti_feproblem.getTimeFromStateArg(state), face.getPoint() + offset);
148 : }
149 : else
150 1686952 : return value(_ti_feproblem.getTimeFromStateArg(state), face.getPoint());
151 : }
152 :
153 : typename Function::ValueType
154 2518466 : Function::evaluate(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
155 : {
156 2518466 : return evaluateHelper(elem_qp, state);
157 : }
158 :
159 : typename Function::ValueType
160 536814 : Function::evaluate(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
161 : {
162 536814 : return evaluateHelper(elem_side_qp, state);
163 : }
164 :
165 : typename Function::ValueType
166 785 : Function::evaluate(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
167 : {
168 785 : return evaluateHelper(elem_point_arg, state);
169 : }
170 :
171 : typename Function::ValueType
172 170713 : Function::evaluate(const NodeArg & node_arg, const Moose::StateArg & state) const
173 : {
174 170713 : return evaluateHelper(node_arg, state);
175 : }
176 :
177 : template <typename R>
178 : typename Function::GradientType
179 7210 : Function::evaluateGradientHelper(const R & r, const Moose::StateArg & state) const
180 : {
181 7210 : return gradient(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
182 : }
183 :
184 : typename Function::GradientType
185 2 : Function::evaluateGradient(const ElemArg & elem_arg, const Moose::StateArg & state) const
186 : {
187 2 : return evaluateGradientHelper(elem_arg, state);
188 : }
189 :
190 : typename Function::GradientType
191 2 : Function::evaluateGradient(const FaceArg & face, const Moose::StateArg & state) const
192 : {
193 2 : return evaluateGradientHelper(face, state);
194 : }
195 :
196 : typename Function::GradientType
197 7202 : Function::evaluateGradient(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
198 : {
199 7202 : return evaluateGradientHelper(elem_qp, state);
200 : }
201 :
202 : typename Function::GradientType
203 2 : Function::evaluateGradient(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
204 : {
205 2 : return evaluateGradientHelper(elem_side_qp, state);
206 : }
207 :
208 : typename Function::GradientType
209 2 : Function::evaluateGradient(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
210 : {
211 2 : return evaluateGradientHelper(elem_point_arg, state);
212 : }
213 :
214 : typename Function::GradientType
215 0 : Function::evaluateGradient(const NodeArg & node_arg, const Moose::StateArg & state) const
216 : {
217 0 : return evaluateGradientHelper(node_arg, state);
218 : }
219 :
220 : template <typename R>
221 : typename Function::DotType
222 2986 : Function::evaluateDotHelper(const R & r, const Moose::StateArg & state) const
223 : {
224 2986 : return timeDerivative(_ti_feproblem.getTimeFromStateArg(state), r.getPoint());
225 : }
226 :
227 : typename Function::DotType
228 290 : Function::evaluateDot(const ElemArg & elem_arg, const Moose::StateArg & state) const
229 : {
230 290 : return evaluateDotHelper(elem_arg, state);
231 : }
232 :
233 : typename Function::DotType
234 2 : Function::evaluateDot(const FaceArg & face, const Moose::StateArg & state) const
235 : {
236 2 : return evaluateDotHelper(face, state);
237 : }
238 :
239 : typename Function::DotType
240 2690 : Function::evaluateDot(const ElemQpArg & elem_qp, const Moose::StateArg & state) const
241 : {
242 2690 : return evaluateDotHelper(elem_qp, state);
243 : }
244 :
245 : typename Function::DotType
246 2 : Function::evaluateDot(const ElemSideQpArg & elem_side_qp, const Moose::StateArg & state) const
247 : {
248 2 : return evaluateDotHelper(elem_side_qp, state);
249 : }
250 :
251 : typename Function::DotType
252 2 : Function::evaluateDot(const ElemPointArg & elem_point_arg, const Moose::StateArg & state) const
253 : {
254 2 : return evaluateDotHelper(elem_point_arg, state);
255 : }
256 :
257 : typename Function::DotType
258 0 : Function::evaluateDot(const NodeArg & node_arg, const Moose::StateArg & state) const
259 : {
260 0 : return evaluateDotHelper(node_arg, state);
261 : }
262 :
263 : void
264 247166 : Function::timestepSetup()
265 : {
266 247166 : FunctorBase<Real>::timestepSetup();
267 247166 : }
268 :
269 : void
270 3058076 : Function::residualSetup()
271 : {
272 3058076 : FunctorBase<Real>::residualSetup();
273 3058076 : }
274 :
275 : void
276 459348 : Function::jacobianSetup()
277 : {
278 459348 : FunctorBase<Real>::jacobianSetup();
279 459348 : }
280 :
281 : void
282 1579199 : Function::customSetup(const ExecFlagType & exec_type)
283 : {
284 1579199 : FunctorBase<Real>::customSetup(exec_type);
285 1579199 : }
|