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 : #pragma once
11 :
12 : #include <tuple>
13 :
14 : #include "MooseFunctorForward.h"
15 : #include "MooseFunctorArguments.h"
16 : #include "FaceArgInterface.h"
17 : #include "MooseMesh.h"
18 : #include "MooseTypes.h"
19 : #include "MooseError.h"
20 : #include "MooseUtils.h"
21 :
22 : #include "libmesh/remote_elem.h"
23 : #include "libmesh/tensor_tools.h"
24 :
25 : #include "metaphysicl/ct_types.h"
26 :
27 : #include <unordered_map>
28 : #include <functional>
29 :
30 : namespace Moose
31 : {
32 : /**
33 : * An enumeration of possible functor evaluation kinds. The available options are value, gradient,
34 : * time derivative (dot), and gradient of time derivative (gradDot)
35 : */
36 : enum class FunctorEvaluationKind
37 : {
38 : Value,
39 : Gradient,
40 : Dot,
41 : GradDot
42 : };
43 :
44 : /**
45 : * A structure that defines the return type of a functor based on the type of the functor and the
46 : * requested evaluation kind, e.g. value, gradient, time derivative, or gradient of time derivative
47 : */
48 : template <typename, FunctorEvaluationKind>
49 : struct FunctorReturnType;
50 :
51 : /**
52 : * The return type for a value evaluation is just the type of the functor
53 : */
54 : template <typename T>
55 : struct FunctorReturnType<T, FunctorEvaluationKind::Value>
56 : {
57 : typedef T type;
58 : };
59 :
60 : /**
61 : * The return type of a gradient evaluation is the rank increment of a value return type. So if the
62 : * value type is Real, then a gradient will be a VectorValue<Real>. This also allows for containers
63 : * of mathematical types. So if a value type is std::vector<Real>, then the gradient type will be
64 : * std::vector<VectorValue<Real>>
65 : */
66 : template <typename T>
67 : struct FunctorReturnType<T, FunctorEvaluationKind::Gradient>
68 : {
69 : typedef typename MetaPhysicL::ReplaceAlgebraicType<
70 : T,
71 : typename libMesh::TensorTools::IncrementRank<
72 : typename MetaPhysicL::ValueType<T>::type>::type>::type type;
73 : };
74 :
75 : /**
76 : * The return type of a time derivative evaluation is the same as the value type
77 : */
78 : template <typename T>
79 : struct FunctorReturnType<T, FunctorEvaluationKind::Dot>
80 : {
81 : typedef T type;
82 : };
83 :
84 : /**
85 : * The return type of a gradient of time derivative evaluation is the same as the gradient type
86 : */
87 : template <typename T>
88 : struct FunctorReturnType<T, FunctorEvaluationKind::GradDot>
89 : {
90 : typedef typename FunctorReturnType<T, FunctorEvaluationKind::Gradient>::type type;
91 : };
92 :
93 : /**
94 : * This structure takes an evaluation kind as a template argument and defines a constant expression
95 : * indicating the associated gradient kind
96 : */
97 : template <FunctorEvaluationKind>
98 : struct FunctorGradientEvaluationKind;
99 :
100 : /**
101 : * The gradient kind associated with a value is simply the gradient
102 : */
103 : template <>
104 : struct FunctorGradientEvaluationKind<FunctorEvaluationKind::Value>
105 : {
106 : static constexpr FunctorEvaluationKind value = FunctorEvaluationKind::Gradient;
107 : };
108 :
109 : /**
110 : * The gradient kind associated with a time derivative is the gradient of the time derivative
111 : */
112 : template <>
113 : struct FunctorGradientEvaluationKind<FunctorEvaluationKind::Dot>
114 : {
115 : static constexpr FunctorEvaluationKind value = FunctorEvaluationKind::GradDot;
116 : };
117 :
118 : /**
119 : * Abstract base class that can be used to hold collections of functors
120 : */
121 : class FunctorAbstract : public FaceArgInterface
122 : {
123 : public:
124 : virtual void residualSetup() = 0;
125 : virtual void jacobianSetup() = 0;
126 : virtual void timestepSetup() = 0;
127 : virtual void customSetup(const ExecFlagType & exec_type) = 0;
128 : };
129 :
130 : /**
131 : * Base class template for functor objects. This class template defines various \p operator()
132 : * overloads that allow a user to evaluate the functor at arbitrary geometric locations. This
133 : * template is meant to enable highly flexible on-the-fly variable and material property
134 : * evaluations
135 : */
136 : template <typename T>
137 : class FunctorBase : public FunctorAbstract
138 : {
139 : public:
140 : using FunctorType = FunctorBase<T>;
141 : using ValueType = T;
142 : /// This rigmarole makes it so that a user can create functors that return containers (std::vector,
143 : /// std::array). This logic will make it such that if a user requests a functor type T that is a
144 : /// container of algebraic types, for example Reals, then the GradientType will be a container of
145 : /// the gradients of those algebraic types, in this example VectorValue<Reals>. So if T is
146 : /// std::vector<Real>, then GradientType will be std::vector<VectorValue<Real>>. As another
147 : /// example: T = std::array<VectorValue<Real>, 1> -> GradientType = std::array<TensorValue<Real>,
148 : /// 1>
149 : using GradientType = typename FunctorReturnType<T, FunctorEvaluationKind::Gradient>::type;
150 : using DotType = ValueType;
151 :
152 1192003 : virtual ~FunctorBase() = default;
153 1242614 : FunctorBase(const MooseFunctorName & name,
154 : const std::set<ExecFlagType> & clearance_schedule = {EXEC_ALWAYS})
155 1242614 : : _always_evaluate(true), _functor_name(name)
156 : {
157 1242614 : setCacheClearanceSchedule(clearance_schedule);
158 1242614 : }
159 :
160 : #ifdef MOOSE_KOKKOS_ENABLED
161 : /**
162 : * Special constructor used for Kokkos functor copy during parallel dispatch
163 : */
164 7379 : FunctorBase(const FunctorBase<T> &, const Moose::Kokkos::FunctorCopy &) {}
165 : #endif
166 :
167 : /**
168 : * Perform a generic evaluation based on the supplied template argument \p FET and supplied
169 : * spatial and temporal arguments
170 : */
171 : template <FunctorEvaluationKind FET, typename Space, typename State>
172 : typename FunctorReturnType<T, FET>::type genericEvaluate(const Space & r,
173 : const State & state) const;
174 :
175 : /// Return the functor name
176 926796 : const MooseFunctorName & functorName() const { return _functor_name; }
177 :
178 : ///@{
179 : /**
180 : * Same as their \p evaluate overloads with the same arguments but allows for caching
181 : * implementation. These are the methods a user will call in their code
182 : */
183 : ValueType operator()(const ElemArg & elem, const StateArg & state) const;
184 : ValueType operator()(const FaceArg & face, const StateArg & state) const;
185 : ValueType operator()(const ElemQpArg & qp, const StateArg & state) const;
186 : ValueType operator()(const ElemSideQpArg & qp, const StateArg & state) const;
187 : ValueType operator()(const ElemPointArg & elem_point, const StateArg & state) const;
188 : ValueType operator()(const NodeArg & node, const StateArg & state) const;
189 : ///@}
190 :
191 : ///@{
192 : /**
193 : * Same as their \p evaluateGradient overloads with the same arguments but allows for caching
194 : * implementation. These are the methods a user will call in their code
195 : */
196 : GradientType gradient(const ElemArg & elem, const StateArg & state) const;
197 : GradientType gradient(const FaceArg & face, const StateArg & state) const;
198 : GradientType gradient(const ElemQpArg & qp, const StateArg & state) const;
199 : GradientType gradient(const ElemSideQpArg & qp, const StateArg & state) const;
200 : GradientType gradient(const ElemPointArg & elem_point, const StateArg & state) const;
201 : GradientType gradient(const NodeArg & node, const StateArg & state) const;
202 : ///@}
203 :
204 : ///@{
205 : /**
206 : * Same as their \p evaluateDot overloads with the same arguments but allows for caching
207 : * implementation. These are the methods a user will call in their code
208 : */
209 : DotType dot(const ElemArg & elem, const StateArg & state) const;
210 : DotType dot(const FaceArg & face, const StateArg & state) const;
211 : DotType dot(const ElemQpArg & qp, const StateArg & state) const;
212 : DotType dot(const ElemSideQpArg & qp, const StateArg & state) const;
213 : DotType dot(const ElemPointArg & elem_point, const StateArg & state) const;
214 : DotType dot(const NodeArg & node, const StateArg & state) const;
215 : ///@}
216 :
217 : ///@{
218 : /**
219 : * Same as their \p evaluateGradDot overloads with the same arguments but allows for caching
220 : * implementation. These are the methods a user will call in their code
221 : */
222 : GradientType gradDot(const ElemArg & elem, const StateArg & state) const;
223 : GradientType gradDot(const FaceArg & face, const StateArg & state) const;
224 : GradientType gradDot(const ElemQpArg & qp, const StateArg & state) const;
225 : GradientType gradDot(const ElemSideQpArg & qp, const StateArg & state) const;
226 : GradientType gradDot(const ElemPointArg & elem_point, const StateArg & state) const;
227 : GradientType gradDot(const NodeArg & node, const StateArg & state) const;
228 : ///@}
229 :
230 : virtual void residualSetup() override;
231 : virtual void jacobianSetup() override;
232 : virtual void timestepSetup() override;
233 : virtual void customSetup(const ExecFlagType & exec_type) override;
234 :
235 : /**
236 : * Set how often to clear the functor evaluation cache
237 : */
238 : void setCacheClearanceSchedule(const std::set<ExecFlagType> & clearance_schedule);
239 :
240 : /**
241 : * Returns whether the functor is defined on this block
242 : */
243 0 : virtual bool hasBlocks(SubdomainID /* id */) const
244 : {
245 0 : mooseError("Block restriction has not been implemented for functor " + functorName());
246 : return false;
247 : }
248 :
249 : /**
250 : * Returns whether this (sided) face is an extrapolated boundary face for
251 : * this functor
252 : */
253 0 : virtual bool isExtrapolatedBoundaryFace(const FaceInfo &, const Elem *, const StateArg &) const
254 : {
255 0 : mooseError("not implemented");
256 : }
257 :
258 : /**
259 : * Returns true if the face is an internal face
260 : */
261 : bool isInternalFace(const FaceInfo &) const;
262 :
263 : /**
264 : * Returns true if this functor is a constant
265 : */
266 0 : virtual bool isConstant() const { return false; }
267 :
268 : virtual bool hasFaceSide(const FaceInfo & fi, const bool fi_elem_side) const override;
269 :
270 : /**
271 : * Examines the incoming face argument. If the face argument producer (residual object,
272 : * postprocessor, etc.) did not indicate a sidedness to the face, e.g. if the \p face_side member
273 : * of the \p FaceArg is \p nullptr, then we may "modify" the sidedness of the argument if we are
274 : * only defined on one side of the face. If the face argument producer \emph has indicated a
275 : * sidedness and we are not defined on that side, then we will error
276 : * @param face The face argument created by the face argument producer, likely a residual object
277 : * @return A face with possibly changed sidedness depending on whether we aren't defined on both
278 : * sides of the face
279 : */
280 : void checkFace(const Moose::FaceArg & face) const;
281 :
282 : /**
283 : * Whether this functor supports evaluation with FaceArg
284 : */
285 : virtual bool supportsFaceArg() const = 0;
286 :
287 : /**
288 : * Whether this functor supports evaluation with ElemSideQpArg
289 : */
290 : virtual bool supportsElemSideQpArg() const = 0;
291 :
292 : protected:
293 : /** @name Functor evaluation routines
294 : * These methods are all for evaluating functors with different kinds of spatial arguments. Each
295 : * of these methods also takes a state argument. For a description of the state argument, please
296 : * see the \p StateArg doxygen
297 : */
298 : ///@{
299 : /**
300 : * Evaluate the functor with a given element. Some example implementations of this method
301 : * could compute an element-average or evaluate at the element centroid
302 : */
303 : virtual ValueType evaluate(const ElemArg & elem, const StateArg & state) const = 0;
304 :
305 : /**
306 : * @param face See the \p FaceArg doxygen
307 : * @param state See the \p StateArg doxygen
308 : * @return The functor evaluated at the requested state and space
309 : */
310 : virtual ValueType evaluate(const FaceArg & face, const StateArg & state) const = 0;
311 :
312 : /**
313 : * @param qp See the \p ElemQpArg doxygen
314 : * @param state See the \p StateArg doxygen
315 : * @return The functor evaluated at the requested state and space
316 : */
317 : virtual ValueType evaluate(const ElemQpArg & qp, const StateArg & state) const = 0;
318 :
319 : /**
320 : * @param side_qp See the \p ElemSideQpArg doxygen
321 : * @param state See the \p StateArg doxygen
322 : * @return The functor evaluated at the requested state and space
323 : */
324 : virtual ValueType evaluate(const ElemSideQpArg & side_qp, const StateArg & state) const = 0;
325 :
326 : /**
327 : * Evaluate the functor with a given element and point. Some example implementations of this
328 : * method could perform a two-term Taylor expansion using cell-centered value and gradient
329 : */
330 : virtual ValueType evaluate(const ElemPointArg & elem_point, const StateArg & state) const = 0;
331 :
332 : virtual ValueType evaluate(const NodeArg & node, const StateArg & state) const = 0;
333 :
334 : /**
335 : * Evaluate the functor gradient with a given element. Some example implementations of this
336 : * method could compute an element-average or evaluate at the element centroid
337 : */
338 2 : virtual GradientType evaluateGradient(const ElemArg &, const StateArg &) const
339 : {
340 2 : mooseError("Element gradient not implemented for functor " + functorName());
341 : }
342 :
343 : /**
344 : * @param face See the \p FaceArg doxygen
345 : * @param state See the \p StateArg doxygen
346 : * @return The functor gradient evaluated at the requested state and space
347 : */
348 2 : virtual GradientType evaluateGradient(const FaceArg &, const StateArg &) const
349 : {
350 2 : mooseError("Face gradient not implemented for functor " + functorName());
351 : }
352 :
353 : /**
354 : * @param qp See the \p ElemQpArg doxygen
355 : * @param state See the \p StateArg doxygen
356 : * @return The functor gradient evaluated at the requested state and space
357 : */
358 2 : virtual GradientType evaluateGradient(const ElemQpArg &, const StateArg &) const
359 : {
360 2 : mooseError("Element quadrature point gradient not implemented for functor " + functorName());
361 : }
362 :
363 : /**
364 : * @param side_qp See the \p ElemSideQpArg doxygen
365 : * @param state See the \p StateArg doxygen
366 : * @return The functor gradient evaluated at the requested state and space
367 : */
368 2 : virtual GradientType evaluateGradient(const ElemSideQpArg &, const StateArg &) const
369 : {
370 2 : mooseError("Element side quadrature point gradient not implemented for functor " +
371 1 : functorName());
372 : }
373 :
374 : /**
375 : * Evaluate the functor gradient with a given element and point
376 : */
377 2 : virtual GradientType evaluateGradient(const ElemPointArg &, const StateArg &) const
378 : {
379 2 : mooseError("Element-point gradient not implemented for functor " + functorName());
380 : }
381 :
382 2 : virtual GradientType evaluateGradient(const NodeArg &, const StateArg &) const
383 : {
384 2 : mooseError("Gradient at node not implemented for functor " + functorName());
385 : }
386 :
387 : /**
388 : * Evaluate the functor time derivative with a given element. Some example implementations of
389 : * this method could compute an element-average or evaluate at the element centroid
390 : */
391 2 : virtual DotType evaluateDot(const ElemArg &, const StateArg &) const
392 : {
393 2 : mooseError("Element time derivative not implemented for functor " + functorName());
394 : }
395 :
396 : /**
397 : * @param face See the \p FaceArg doxygen
398 : * @param state See the \p StateArg doxygen
399 : * @return The functor time derivative evaluated at the requested state and space
400 : */
401 2 : virtual DotType evaluateDot(const FaceArg &, const StateArg &) const
402 : {
403 2 : mooseError("Face time derivative not implemented for functor " + functorName());
404 : }
405 :
406 : /**
407 : * @param qp See the \p ElemQpArg doxygen
408 : * @param state See the \p StateArg doxygen
409 : * @return The functor time derivative evaluated at the requested state and space
410 : */
411 6 : virtual DotType evaluateDot(const ElemQpArg &, const StateArg &) const
412 : {
413 6 : mooseError("Element quadrature point time derivative not implemented for functor " +
414 5 : functorName());
415 : }
416 :
417 : /**
418 : * @param side_qp See the \p ElemSideQpArg doxygen
419 : * @param state See the \p StateArg doxygen
420 : * @return The functor time derivative evaluated at the requested state and space
421 : */
422 2 : virtual DotType evaluateDot(const ElemSideQpArg &, const StateArg &) const
423 : {
424 2 : mooseError("Element side quadrature point time derivative not implemented for functor " +
425 1 : functorName());
426 : }
427 :
428 : /**
429 : * Evaluate the functor time derivative with a given element and point
430 : */
431 2 : virtual DotType evaluateDot(const ElemPointArg &, const StateArg &) const
432 : {
433 2 : mooseError("Element-point time derivative not implemented for functor " + functorName());
434 : }
435 :
436 2 : virtual DotType evaluateDot(const NodeArg &, const StateArg &) const
437 : {
438 2 : mooseError("Time derivative at node not implemented for functor " + functorName());
439 : }
440 :
441 : /**
442 : * Evaluate the functor gradient-dot with a given element. Some example implementations of this
443 : * method could compute an element-average or evaluate at the element centroid
444 : */
445 0 : virtual GradientType evaluateGradDot(const ElemArg &, const StateArg &) const
446 : {
447 0 : mooseError("Element gradient-dot not implemented for functor " + functorName());
448 : }
449 :
450 : /**
451 : * @param face See the \p FaceArg doxygen
452 : * @param state See the \p StateArg doxygen
453 : * @return The functor gradient-dot evaluated at the requested state and space
454 : */
455 0 : virtual GradientType evaluateGradDot(const FaceArg &, const StateArg &) const
456 : {
457 0 : mooseError("Face gradient-dot not implemented for functor " + functorName());
458 : }
459 :
460 : /**
461 : * @param qp See the \p ElemQpArg doxygen
462 : * @param state See the \p StateArg doxygen
463 : * @return The functor gradient-dot evaluated at the requested state and space
464 : */
465 0 : virtual GradientType evaluateGradDot(const ElemQpArg &, const StateArg &) const
466 : {
467 0 : mooseError("Element quadrature point gradient-dot not implemented for functor " +
468 0 : functorName());
469 : }
470 :
471 : /**
472 : * @param side_qp See the \p ElemSideQpArg doxygen
473 : * @param state See the \p StateArg doxygen
474 : * @return The functor gradient-dot evaluated at the requested state and space
475 : */
476 0 : virtual GradientType evaluateGradDot(const ElemSideQpArg &, const StateArg &) const
477 : {
478 0 : mooseError("Element side quadrature point gradient-dot not implemented for functor " +
479 0 : functorName());
480 : }
481 :
482 : /**
483 : * Evaluate the functor gradient-dot with a given element and point
484 : */
485 0 : virtual GradientType evaluateGradDot(const ElemPointArg &, const StateArg &) const
486 : {
487 0 : mooseError("Element-point gradient-dot not implemented for functor " + functorName());
488 : }
489 :
490 0 : virtual GradientType evaluateGradDot(const NodeArg &, const StateArg &) const
491 : {
492 0 : mooseError("Gradient-dot at node not implemented for functor " + functorName());
493 : }
494 : ///@}
495 :
496 : private:
497 : /**
498 : * clear cache data
499 : */
500 : void clearCacheData();
501 :
502 : /**
503 : * check a qp cache and if invalid then evaluate
504 : */
505 : template <typename SpaceArg, typename StateArg>
506 : ValueType queryQpCache(unsigned int qp,
507 : const libMesh::QBase & qrule,
508 : std::vector<std::pair<bool, T>> & qp_cache_data,
509 : const SpaceArg & space,
510 : const StateArg & state) const;
511 :
512 : /**
513 : * check a finite volume spatial argument cache and if invalid then evaluate
514 : */
515 : template <typename SpaceArg>
516 : ValueType queryFVArgCache(std::map<SpaceArg, ValueType> & cache_data,
517 : const SpaceArg & space) const;
518 :
519 : /// How often to clear the material property cache
520 : std::set<ExecFlagType> _clearance_schedule;
521 :
522 : /// Boolean to check if we always need evaluation
523 : bool _always_evaluate;
524 :
525 : // Data for traditional element-quadrature point property evaluations which are useful for
526 : // caching implementation
527 :
528 : /// Current key for qp map cache
529 : mutable dof_id_type _current_qp_map_key = libMesh::DofObject::invalid_id;
530 :
531 : /// Current value for qp map cache
532 : mutable std::vector<std::pair<bool, ValueType>> * _current_qp_map_value = nullptr;
533 :
534 : /// Cached element quadrature point functor property evaluations. The map key is the element
535 : /// id. The map values should have size corresponding to the number of quadrature points on the
536 : /// element. The vector elements are pairs. The first member of the pair indicates whether a
537 : /// cached value has been computed. The second member of the pair is the (cached) value. If the
538 : /// boolean is false, then the value cannot be trusted
539 : mutable std::unordered_map<dof_id_type, std::vector<std::pair<bool, ValueType>>> _qp_to_value;
540 :
541 : // Data for traditional element-side-quadrature point property evaluations which are useful for
542 : // caching implementation
543 :
544 : /// Current key for side-qp map cache
545 : mutable dof_id_type _current_side_qp_map_key = libMesh::DofObject::invalid_id;
546 :
547 : /// Current value for side-qp map cache
548 : mutable std::vector<std::vector<std::pair<bool, ValueType>>> * _current_side_qp_map_value =
549 : nullptr;
550 :
551 : /// Cached element quadrature point functor property evaluations. The map key is the element
552 : /// id. The map values are a multi-dimensional vector (or vector of vectors) with the first index
553 : /// corresponding to the side and the second index corresponding to the quadrature point
554 : /// index. The elements returned after double indexing are pairs. The first member of the pair
555 : /// indicates whether a cached value has been computed. The second member of the pair is the
556 : /// (cached) value. If the boolean is false, then the value cannot be trusted
557 : mutable std::unordered_map<dof_id_type, std::vector<std::vector<std::pair<bool, ValueType>>>>
558 : _side_qp_to_value;
559 :
560 : /// Map from element arguments to their cached evaluations
561 : mutable std::map<ElemArg, ValueType> _elem_arg_to_value;
562 :
563 : /// Map from face arguments to their cached evaluations
564 : mutable std::map<FaceArg, ValueType> _face_arg_to_value;
565 :
566 : /// Map from nodal arguments to their cached evaluations
567 : mutable std::map<NodeArg, ValueType> _node_arg_to_value;
568 :
569 : /// name of the functor
570 : MooseFunctorName _functor_name;
571 : };
572 :
573 : template <typename T>
574 : bool
575 175238108 : FunctorBase<T>::isInternalFace(const FaceInfo & fi) const
576 : {
577 175238108 : if (!fi.neighborPtr())
578 3594014 : return false;
579 :
580 171644094 : return hasBlocks(fi.elem().subdomain_id()) && hasBlocks(fi.neighborPtr()->subdomain_id());
581 : }
582 :
583 : template <typename T>
584 : template <typename SpaceArg>
585 : typename FunctorBase<T>::ValueType
586 0 : FunctorBase<T>::queryFVArgCache(std::map<SpaceArg, ValueType> & cache_data,
587 : const SpaceArg & space) const
588 : {
589 : // We don't want to evaluate if the key already exists, so instead we value initialize
590 0 : auto [it, inserted] = cache_data.try_emplace(space, ValueType());
591 0 : auto & value = it->second;
592 :
593 0 : if (inserted)
594 : // value not ready to go
595 : // this function is only called from functions that assert we are in the current time state
596 0 : value = evaluate(space, currentState());
597 :
598 0 : return value;
599 : }
600 :
601 : template <typename T>
602 : typename FunctorBase<T>::ValueType
603 228177949 : FunctorBase<T>::operator()(const ElemArg & elem, const StateArg & state) const
604 : {
605 228177949 : if (_always_evaluate)
606 228177949 : return evaluate(elem, state);
607 :
608 : mooseAssert(state.state == 0,
609 : "Cached evaluations are only currently supported for the current state.");
610 :
611 0 : return queryFVArgCache(_elem_arg_to_value, elem);
612 : }
613 :
614 : template <typename T>
615 : typename FunctorBase<T>::ValueType
616 79864891 : FunctorBase<T>::operator()(const FaceArg & face_in, const StateArg & state) const
617 : {
618 79864891 : checkFace(face_in);
619 :
620 79864891 : if (_always_evaluate)
621 79864891 : return evaluate(face_in, state);
622 :
623 : mooseAssert(state.state == 0,
624 : "Cached evaluations are only currently supported for the current state.");
625 :
626 0 : return queryFVArgCache(_face_arg_to_value, face_in);
627 : }
628 :
629 : template <typename T>
630 : template <typename SpaceArg, typename StateArg>
631 : typename FunctorBase<T>::ValueType
632 0 : FunctorBase<T>::queryQpCache(const unsigned int qp,
633 : const libMesh::QBase & qrule,
634 : std::vector<std::pair<bool, ValueType>> & qp_cache_data,
635 : const SpaceArg & space,
636 : const StateArg & state) const
637 : {
638 : // Check and see whether we even have sized for this quadrature point. If we haven't then we
639 : // must evaluate
640 0 : if (qp >= qp_cache_data.size())
641 : {
642 0 : qp_cache_data.resize(qrule.n_points(), std::make_pair(false, ValueType()));
643 0 : auto & pr = qp_cache_data[qp];
644 0 : pr.second = evaluate(space, state);
645 0 : pr.first = true;
646 0 : return pr.second;
647 : }
648 :
649 : // We've already sized for this qp, so let's see whether we have a valid cache value
650 0 : auto & pr = qp_cache_data[qp];
651 0 : if (pr.first)
652 0 : return pr.second;
653 :
654 : // No valid cache value so evaluate
655 0 : pr.second = evaluate(space, state);
656 0 : pr.first = true;
657 0 : return pr.second;
658 : }
659 :
660 : template <typename T>
661 : typename FunctorBase<T>::ValueType
662 31613570 : FunctorBase<T>::operator()(const ElemQpArg & elem_qp, const StateArg & state) const
663 : {
664 31613570 : if (_always_evaluate)
665 31613570 : return evaluate(elem_qp, state);
666 :
667 0 : const auto elem_id = elem_qp.elem->id();
668 0 : if (elem_id != _current_qp_map_key)
669 : {
670 0 : _current_qp_map_key = elem_id;
671 0 : _current_qp_map_value = &_qp_to_value[elem_id];
672 : }
673 0 : auto & qp_data = *_current_qp_map_value;
674 0 : const auto qp = elem_qp.qp;
675 0 : const auto * const qrule = elem_qp.qrule;
676 : mooseAssert(qrule, "qrule must be non-null");
677 :
678 0 : return queryQpCache(qp, *qrule, qp_data, elem_qp, state);
679 : }
680 :
681 : template <typename T>
682 : typename FunctorBase<T>::ValueType
683 1273268 : FunctorBase<T>::operator()(const ElemSideQpArg & elem_side_qp, const StateArg & state) const
684 : {
685 1273268 : if (_always_evaluate)
686 1273268 : return evaluate(elem_side_qp, state);
687 :
688 0 : const Elem * const elem = elem_side_qp.elem;
689 : mooseAssert(elem, "elem must be non-null");
690 0 : const auto elem_id = elem->id();
691 0 : if (elem_id != _current_side_qp_map_key)
692 : {
693 0 : _current_side_qp_map_key = elem_id;
694 0 : _current_side_qp_map_value = &_side_qp_to_value[elem_id];
695 : }
696 0 : auto & side_qp_data = *_current_side_qp_map_value;
697 0 : const auto side = elem_side_qp.side;
698 0 : const auto qp = elem_side_qp.qp;
699 0 : const auto * const qrule = elem_side_qp.qrule;
700 : mooseAssert(qrule, "qrule must be non-null");
701 :
702 : // Check and see whether we even have sized for this side
703 0 : if (side >= side_qp_data.size())
704 0 : side_qp_data.resize(elem->n_sides());
705 :
706 : // Ok we were sized enough for our side
707 0 : auto & qp_data = side_qp_data[side];
708 0 : return queryQpCache(qp, *qrule, qp_data, elem_side_qp, state);
709 : }
710 :
711 : template <typename T>
712 : typename FunctorBase<T>::ValueType
713 183444 : FunctorBase<T>::operator()(const ElemPointArg & elem_point, const StateArg & state) const
714 : {
715 183444 : return evaluate(elem_point, state);
716 : }
717 :
718 : template <typename T>
719 : void
720 1242614 : FunctorBase<T>::setCacheClearanceSchedule(const std::set<ExecFlagType> & clearance_schedule)
721 : {
722 1242614 : if (clearance_schedule.count(EXEC_ALWAYS))
723 1242140 : _always_evaluate = true;
724 :
725 1242614 : _clearance_schedule = clearance_schedule;
726 1242614 : }
727 :
728 : template <typename T>
729 : typename FunctorBase<T>::ValueType
730 790113 : FunctorBase<T>::operator()(const NodeArg & node, const StateArg & state) const
731 : {
732 : mooseAssert(node.subdomain_ids, "Subdomain IDs must be supplied to the node argument");
733 790113 : return evaluate(node, state);
734 : }
735 :
736 : template <typename T>
737 : void
738 84573984 : FunctorBase<T>::checkFace(const Moose::FaceArg &
739 : #if DEBUG
740 : face
741 : #endif
742 : ) const
743 : {
744 : #if DEBUG
745 : const Elem * const elem = face.face_side;
746 : const FaceInfo * const fi = face.fi;
747 : mooseAssert(fi, "face info should be non-null");
748 : bool check_elem_def = false;
749 : bool check_neighbor_def = false;
750 : // We check if the functor is defined on both sides of the face
751 : if (!elem)
752 : {
753 : if (!hasFaceSide(*fi, true))
754 : check_neighbor_def = true;
755 : else if (!hasFaceSide(*fi, false))
756 : check_elem_def = true;
757 : }
758 : else if (elem == fi->elemPtr())
759 : check_elem_def = true;
760 : else
761 : {
762 : mooseAssert(elem == fi->neighborPtr(), "This has to match something");
763 : check_neighbor_def = true;
764 : }
765 :
766 : if (check_elem_def && !hasFaceSide(*fi, true))
767 : {
768 : std::string additional_message = "It is not defined on the neighbor side either.";
769 : if (hasFaceSide(*fi, false))
770 : additional_message = "It is however defined on the neighbor side.";
771 : additional_message += " Face centroid: " + Moose::stringify(fi->faceCentroid());
772 : mooseError(_functor_name,
773 : " is not defined on the element side of the face information, but a face argument "
774 : "producer "
775 : "(e.g. residual object, postprocessor, etc.) has requested evaluation there.\n",
776 : additional_message);
777 : }
778 : if (check_neighbor_def && !hasFaceSide(*fi, false))
779 : {
780 : std::string additional_message = "It is not defined on the element side either.";
781 : if (hasFaceSide(*fi, true))
782 : additional_message = "It is however defined on the element side.";
783 : additional_message += " Face centroid: " + Moose::stringify(fi->faceCentroid());
784 : mooseError(
785 : _functor_name,
786 : " is not defined on the neighbor side of the face information, but a face argument "
787 : "producer (e.g. residual object, postprocessor, etc.) has requested evaluation there.\n",
788 : additional_message);
789 : }
790 : #endif
791 84573984 : }
792 :
793 : template <typename T>
794 : void
795 29270 : FunctorBase<T>::clearCacheData()
796 : {
797 29270 : for (auto & map_pr : _qp_to_value)
798 0 : for (auto & pr : map_pr.second)
799 0 : pr.first = false;
800 :
801 29270 : for (auto & map_pr : _side_qp_to_value)
802 : {
803 0 : auto & side_vector = map_pr.second;
804 0 : for (auto & qp_vector : side_vector)
805 0 : for (auto & pr : qp_vector)
806 0 : pr.first = false;
807 : }
808 :
809 29270 : _current_qp_map_key = libMesh::DofObject::invalid_id;
810 29270 : _current_qp_map_value = nullptr;
811 29270 : _current_side_qp_map_key = libMesh::DofObject::invalid_id;
812 29270 : _current_side_qp_map_value = nullptr;
813 :
814 29270 : _elem_arg_to_value.clear();
815 29270 : _face_arg_to_value.clear();
816 29270 : _node_arg_to_value.clear();
817 29270 : }
818 :
819 : template <typename T>
820 : void
821 1399419 : FunctorBase<T>::timestepSetup()
822 : {
823 1399419 : if (_clearance_schedule.count(EXEC_TIMESTEP_BEGIN))
824 25 : clearCacheData();
825 1399419 : }
826 :
827 : template <typename T>
828 : void
829 10252831 : FunctorBase<T>::residualSetup()
830 : {
831 10252831 : if (_clearance_schedule.count(EXEC_LINEAR))
832 28063 : clearCacheData();
833 10252831 : }
834 :
835 : template <typename T>
836 : void
837 1658494 : FunctorBase<T>::jacobianSetup()
838 : {
839 1658494 : if (_clearance_schedule.count(EXEC_NONLINEAR))
840 1182 : clearCacheData();
841 1658494 : }
842 :
843 : template <typename T>
844 : void
845 1515716 : FunctorBase<T>::customSetup(const ExecFlagType & exec_type)
846 : {
847 1515716 : if (_clearance_schedule.count(exec_type))
848 0 : clearCacheData();
849 1515716 : }
850 :
851 : template <typename T>
852 : typename FunctorBase<T>::GradientType
853 4960783 : FunctorBase<T>::gradient(const ElemArg & elem, const StateArg & state) const
854 : {
855 4960783 : return evaluateGradient(elem, state);
856 : }
857 :
858 : template <typename T>
859 : typename FunctorBase<T>::GradientType
860 4709051 : FunctorBase<T>::gradient(const FaceArg & face, const StateArg & state) const
861 : {
862 4709051 : checkFace(face);
863 4709051 : return evaluateGradient(face, state);
864 : }
865 :
866 : template <typename T>
867 : typename FunctorBase<T>::GradientType
868 11499482 : FunctorBase<T>::gradient(const ElemQpArg & elem_qp, const StateArg & state) const
869 : {
870 11499482 : return evaluateGradient(elem_qp, state);
871 : }
872 :
873 : template <typename T>
874 : typename FunctorBase<T>::GradientType
875 22 : FunctorBase<T>::gradient(const ElemSideQpArg & elem_side_qp, const StateArg & state) const
876 : {
877 22 : return evaluateGradient(elem_side_qp, state);
878 : }
879 :
880 : template <typename T>
881 : typename FunctorBase<T>::GradientType
882 22 : FunctorBase<T>::gradient(const ElemPointArg & elem_point, const StateArg & state) const
883 : {
884 22 : return evaluateGradient(elem_point, state);
885 : }
886 :
887 : template <typename T>
888 : typename FunctorBase<T>::GradientType
889 50 : FunctorBase<T>::gradient(const NodeArg & node, const StateArg & state) const
890 : {
891 50 : return evaluateGradient(node, state);
892 : }
893 :
894 : template <typename T>
895 : typename FunctorBase<T>::DotType
896 76314 : FunctorBase<T>::dot(const ElemArg & elem, const StateArg & state) const
897 : {
898 76314 : return evaluateDot(elem, state);
899 : }
900 :
901 : template <typename T>
902 : typename FunctorBase<T>::DotType
903 22 : FunctorBase<T>::dot(const FaceArg & face, const StateArg & state) const
904 : {
905 22 : checkFace(face);
906 22 : return evaluateDot(face, state);
907 : }
908 :
909 : template <typename T>
910 : typename FunctorBase<T>::DotType
911 843314 : FunctorBase<T>::dot(const ElemQpArg & elem_qp, const StateArg & state) const
912 : {
913 843314 : return evaluateDot(elem_qp, state);
914 : }
915 :
916 : template <typename T>
917 : typename FunctorBase<T>::DotType
918 22 : FunctorBase<T>::dot(const ElemSideQpArg & elem_side_qp, const StateArg & state) const
919 : {
920 22 : return evaluateDot(elem_side_qp, state);
921 : }
922 :
923 : template <typename T>
924 : typename FunctorBase<T>::DotType
925 22 : FunctorBase<T>::dot(const ElemPointArg & elem_point, const StateArg & state) const
926 : {
927 22 : return evaluateDot(elem_point, state);
928 : }
929 :
930 : template <typename T>
931 : typename FunctorBase<T>::DotType
932 18 : FunctorBase<T>::dot(const NodeArg & node, const StateArg & state) const
933 : {
934 18 : return evaluateDot(node, state);
935 : }
936 :
937 : template <typename T>
938 : typename FunctorBase<T>::GradientType
939 2184 : FunctorBase<T>::gradDot(const ElemArg & elem, const StateArg & state) const
940 : {
941 2184 : return evaluateGradDot(elem, state);
942 : }
943 :
944 : template <typename T>
945 : typename FunctorBase<T>::GradientType
946 20 : FunctorBase<T>::gradDot(const FaceArg & face, const StateArg & state) const
947 : {
948 20 : checkFace(face);
949 20 : return evaluateGradDot(face, state);
950 : }
951 :
952 : template <typename T>
953 : typename FunctorBase<T>::GradientType
954 16 : FunctorBase<T>::gradDot(const ElemQpArg & elem_qp, const StateArg & state) const
955 : {
956 16 : return evaluateGradDot(elem_qp, state);
957 : }
958 :
959 : template <typename T>
960 : typename FunctorBase<T>::GradientType
961 16 : FunctorBase<T>::gradDot(const ElemSideQpArg & elem_side_qp, const StateArg & state) const
962 : {
963 16 : return evaluateGradDot(elem_side_qp, state);
964 : }
965 :
966 : template <typename T>
967 : typename FunctorBase<T>::GradientType
968 16 : FunctorBase<T>::gradDot(const ElemPointArg & elem_point, const StateArg & state) const
969 : {
970 16 : return evaluateGradDot(elem_point, state);
971 : }
972 :
973 : template <typename T>
974 : typename FunctorBase<T>::GradientType
975 16 : FunctorBase<T>::gradDot(const NodeArg & node, const StateArg & state) const
976 : {
977 16 : return evaluateGradDot(node, state);
978 : }
979 :
980 : template <typename T>
981 : bool
982 4217122 : FunctorBase<T>::hasFaceSide(const FaceInfo & fi, const bool fi_elem_side) const
983 : {
984 4217122 : if (fi_elem_side)
985 2108645 : return hasBlocks(fi.elem().subdomain_id());
986 : else
987 2108477 : return fi.neighborPtr() && hasBlocks(fi.neighbor().subdomain_id());
988 : }
989 :
990 : template <typename T>
991 : template <FunctorEvaluationKind FET, typename Space, typename State>
992 : typename FunctorReturnType<T, FET>::type
993 166944438 : FunctorBase<T>::genericEvaluate(const Space & r, const State & state) const
994 : {
995 : if constexpr (FET == FunctorEvaluationKind::Value)
996 157527596 : return (*this)(r, state);
997 : else if constexpr (FET == FunctorEvaluationKind::Gradient)
998 9416842 : return gradient(r, state);
999 : else if constexpr (FET == FunctorEvaluationKind::Dot)
1000 0 : return dot(r, state);
1001 : else
1002 0 : return gradDot(r, state);
1003 : }
1004 :
1005 : /**
1006 : * A non-templated base class for functors that allow an owner object to hold
1007 : * different class template instantiations of \p Functor in a single container
1008 : */
1009 : class FunctorEnvelopeBase
1010 : {
1011 : public:
1012 621194 : FunctorEnvelopeBase() = default;
1013 592138 : virtual ~FunctorEnvelopeBase() = default;
1014 :
1015 : /**
1016 : * @return Whether this envelope wraps a null functor
1017 : */
1018 : virtual bool wrapsNull() const = 0;
1019 :
1020 : /**
1021 : * @return The return type, as a string, of the functor this envelope wraps
1022 : */
1023 : virtual std::string returnType() const = 0;
1024 :
1025 : /**
1026 : * @return Whether this envelope wraps a constant functor
1027 : */
1028 : virtual bool isConstant() const = 0;
1029 :
1030 : /**
1031 : * @return Whether this envelope owns its wrapped functor. This envelope may briefly own null
1032 : * functors during simulation setup or it may own non-AD or AD wrappers of "true" functors, but we
1033 : * should never own any "true" functors, e.g. we expect memory of "true" functors to be managed by
1034 : * warehouses (e.g. variable, function, etc.), or by the \p SubProblem itself. With this
1035 : * expectation, we don't have to worry about performing setup calls
1036 : */
1037 : virtual bool ownsWrappedFunctor() const = 0;
1038 : };
1039 :
1040 : /**
1041 : * This is a wrapper that forwards calls to the implementation,
1042 : * which can be switched out at any time without disturbing references to
1043 : * FunctorBase. Implementation motivated by https://stackoverflow.com/a/65455485/4493669
1044 : */
1045 : template <typename T>
1046 : class FunctorEnvelope final : public FunctorBase<T>, public FunctorEnvelopeBase
1047 : {
1048 : public:
1049 : using typename Moose::FunctorBase<T>::ValueType;
1050 : using typename Moose::FunctorBase<T>::GradientType;
1051 : using typename Moose::FunctorBase<T>::DotType;
1052 :
1053 : /**
1054 : * @param wrapped The functor to wrap. We will *not* not own the wrapped object
1055 : */
1056 304946 : FunctorEnvelope(const FunctorBase<T> & wrapped)
1057 914838 : : FunctorBase<T>("wraps_" + wrapped.functorName()), FunctorEnvelopeBase(), _wrapped(&wrapped)
1058 : {
1059 609892 : }
1060 :
1061 : /**
1062 : * @param wrapped A unique pointer around the functor to wrap. We *will* own the wrapped object,
1063 : * e.g. if we are ever destructed or we are reassigned to wrap another functor, then this functor
1064 : * will be destructed
1065 : */
1066 316248 : FunctorEnvelope(std::unique_ptr<FunctorBase<T>> && wrapped)
1067 316248 : : FunctorBase<T>("wraps_" + wrapped->functorName()),
1068 : FunctorEnvelopeBase(),
1069 316248 : _owned(std::move(wrapped)),
1070 1581240 : _wrapped(_owned.get())
1071 : {
1072 632496 : }
1073 :
1074 : /**
1075 : * Prevent wrapping of a temporary object. If we are to own a functor, the unique_ptr constructor
1076 : * overload should be used
1077 : */
1078 : FunctorEnvelope(FunctorBase<T> &&) = delete;
1079 :
1080 : /**
1081 : * @param wrapped The functor to wrap. We will *not* not own the wrapped object. If we previously
1082 : * owned a functor, it will be destructed
1083 : */
1084 135 : void assign(const FunctorBase<T> & wrapped)
1085 : {
1086 135 : _owned.reset();
1087 135 : _wrapped = &wrapped;
1088 135 : }
1089 :
1090 : /**
1091 : * @param wrapped A unique pointer around the functor to wrap. We *will* own the wrapped object.
1092 : * If we previously owned a functor, it will be destructed
1093 : */
1094 135 : void assign(std::unique_ptr<FunctorBase<T>> && wrapped)
1095 : {
1096 135 : _owned = std::move(wrapped);
1097 135 : _wrapped = _owned.get();
1098 135 : }
1099 :
1100 : /**
1101 : * Prevent wrapping of a temporary object. If we are to own a functor, the unique_ptr assign
1102 : * overload should be used
1103 : */
1104 : void assign(FunctorBase<T> &&) = delete;
1105 :
1106 : FunctorEnvelope(const FunctorEnvelope &) = delete;
1107 : FunctorEnvelope(FunctorEnvelope &&) = delete;
1108 : FunctorEnvelope & operator=(const FunctorEnvelope &) = delete;
1109 : FunctorEnvelope & operator=(FunctorEnvelope &&) = delete;
1110 :
1111 1184274 : virtual ~FunctorEnvelope() = default;
1112 :
1113 : /**
1114 : * @return whether this object wraps a null functor
1115 : */
1116 297794 : virtual bool wrapsNull() const override { return wrapsType<NullFunctor<T>>(); }
1117 :
1118 : /**
1119 : * @return a string representation of the return type of this functor
1120 : */
1121 0 : virtual std::string returnType() const override { return libMesh::demangle(typeid(T).name()); }
1122 :
1123 297794 : virtual bool ownsWrappedFunctor() const override { return _owned.get(); }
1124 :
1125 : /**
1126 : * @return whether the wrapped object is of the requested type
1127 : */
1128 : template <typename T2>
1129 323197 : bool wrapsType() const
1130 : {
1131 323197 : return dynamic_cast<const T2 *>(_wrapped);
1132 : }
1133 :
1134 0 : virtual bool isExtrapolatedBoundaryFace(const FaceInfo & fi,
1135 : const Elem * const elem,
1136 : const StateArg & state) const override
1137 : {
1138 0 : return _wrapped->isExtrapolatedBoundaryFace(fi, elem, state);
1139 : }
1140 72 : virtual bool isConstant() const override { return _wrapped->isConstant(); }
1141 547578 : virtual bool hasBlocks(const SubdomainID id) const override { return _wrapped->hasBlocks(id); }
1142 4211780 : virtual bool hasFaceSide(const FaceInfo & fi, const bool fi_elem_side) const override
1143 : {
1144 4211780 : return _wrapped->hasFaceSide(fi, fi_elem_side);
1145 : }
1146 :
1147 1376 : bool supportsFaceArg() const override final { return true; }
1148 588 : bool supportsElemSideQpArg() const override final { return true; }
1149 :
1150 : protected:
1151 : ///@{
1152 : /**
1153 : * Forward calls to wrapped object
1154 : */
1155 25001830 : virtual ValueType evaluate(const ElemArg & elem, const StateArg & state) const override
1156 : {
1157 25001830 : return _wrapped->operator()(elem, state);
1158 : }
1159 3723113 : virtual ValueType evaluate(const FaceArg & face, const StateArg & state) const override
1160 : {
1161 3723113 : return _wrapped->operator()(face, state);
1162 : }
1163 15683110 : virtual ValueType evaluate(const ElemQpArg & qp, const StateArg & state) const override
1164 : {
1165 15683110 : return _wrapped->operator()(qp, state);
1166 : }
1167 631470 : virtual ValueType evaluate(const ElemSideQpArg & qp, const StateArg & state) const override
1168 : {
1169 631470 : return _wrapped->operator()(qp, state);
1170 : }
1171 75419 : virtual ValueType evaluate(const ElemPointArg & elem_point, const StateArg & state) const override
1172 : {
1173 75419 : return _wrapped->operator()(elem_point, state);
1174 : }
1175 278614 : virtual ValueType evaluate(const NodeArg & node, const StateArg & state) const override
1176 : {
1177 278614 : return _wrapped->operator()(node, state);
1178 : }
1179 :
1180 32378 : virtual GradientType evaluateGradient(const ElemArg & elem, const StateArg & state) const override
1181 : {
1182 32378 : return _wrapped->gradient(elem, state);
1183 : }
1184 2 : virtual GradientType evaluateGradient(const FaceArg & face, const StateArg & state) const override
1185 : {
1186 2 : return _wrapped->gradient(face, state);
1187 : }
1188 75130 : virtual GradientType evaluateGradient(const ElemQpArg & qp, const StateArg & state) const override
1189 : {
1190 75130 : return _wrapped->gradient(qp, state);
1191 : }
1192 2 : virtual GradientType evaluateGradient(const ElemSideQpArg & qp,
1193 : const StateArg & state) const override
1194 : {
1195 2 : return _wrapped->gradient(qp, state);
1196 : }
1197 2 : virtual GradientType evaluateGradient(const ElemPointArg & elem_point,
1198 : const StateArg & state) const override
1199 : {
1200 2 : return _wrapped->gradient(elem_point, state);
1201 : }
1202 16 : virtual GradientType evaluateGradient(const NodeArg & node, const StateArg & state) const override
1203 : {
1204 16 : return _wrapped->gradient(node, state);
1205 : }
1206 :
1207 9874 : virtual DotType evaluateDot(const ElemArg & elem, const StateArg & state) const override
1208 : {
1209 9874 : return _wrapped->dot(elem, state);
1210 : }
1211 2 : virtual DotType evaluateDot(const FaceArg & face, const StateArg & state) const override
1212 : {
1213 2 : return _wrapped->dot(face, state);
1214 : }
1215 2844 : virtual DotType evaluateDot(const ElemQpArg & qp, const StateArg & state) const override
1216 : {
1217 2844 : return _wrapped->dot(qp, state);
1218 : }
1219 2 : virtual DotType evaluateDot(const ElemSideQpArg & qp, const StateArg & state) const override
1220 : {
1221 2 : return _wrapped->dot(qp, state);
1222 : }
1223 2 : virtual DotType evaluateDot(const ElemPointArg & elem_point,
1224 : const StateArg & state) const override
1225 : {
1226 2 : return _wrapped->dot(elem_point, state);
1227 : }
1228 0 : virtual DotType evaluateDot(const NodeArg & node, const StateArg & state) const override
1229 : {
1230 0 : return _wrapped->dot(node, state);
1231 : }
1232 :
1233 1080 : virtual GradientType evaluateGradDot(const ElemArg & elem, const StateArg & state) const override
1234 : {
1235 1080 : return _wrapped->gradDot(elem, state);
1236 : }
1237 0 : virtual GradientType evaluateGradDot(const FaceArg & face, const StateArg & state) const override
1238 : {
1239 0 : return _wrapped->gradDot(face, state);
1240 : }
1241 0 : virtual GradientType evaluateGradDot(const ElemQpArg & qp, const StateArg & state) const override
1242 : {
1243 0 : return _wrapped->gradDot(qp, state);
1244 : }
1245 0 : virtual GradientType evaluateGradDot(const ElemSideQpArg & qp,
1246 : const StateArg & state) const override
1247 : {
1248 0 : return _wrapped->gradDot(qp, state);
1249 : }
1250 0 : virtual GradientType evaluateGradDot(const ElemPointArg & elem_point,
1251 : const StateArg & state) const override
1252 : {
1253 0 : return _wrapped->gradDot(elem_point, state);
1254 : }
1255 0 : virtual GradientType evaluateGradDot(const NodeArg & node, const StateArg & state) const override
1256 : {
1257 0 : return _wrapped->gradDot(node, state);
1258 : }
1259 : ///@}
1260 :
1261 : private:
1262 : /// Our wrapped object
1263 : std::unique_ptr<FunctorBase<T>> _owned;
1264 : const FunctorBase<T> * _wrapped;
1265 :
1266 : friend class ::SubProblem;
1267 : };
1268 :
1269 : /**
1270 : * Class template for creating constant functors
1271 : */
1272 : template <typename T>
1273 : class ConstantFunctor final : public FunctorBase<T>
1274 : {
1275 : public:
1276 : using typename FunctorBase<T>::FunctorType;
1277 : using typename FunctorBase<T>::ValueType;
1278 : using typename FunctorBase<T>::GradientType;
1279 : using typename FunctorBase<T>::DotType;
1280 :
1281 4706 : ConstantFunctor(const ValueType & value)
1282 14118 : : FunctorBase<T>("constant_" + std::to_string(value)), _value(value)
1283 : {
1284 9412 : }
1285 6330 : ConstantFunctor(ValueType && value)
1286 18990 : : FunctorBase<T>("constant_" + std::to_string(MetaPhysicL::raw_value(value))), _value(value)
1287 : {
1288 12660 : }
1289 :
1290 72 : virtual bool isConstant() const override { return true; }
1291 :
1292 1677308 : bool hasBlocks(SubdomainID /* id */) const override { return true; }
1293 :
1294 0 : bool supportsFaceArg() const override final { return true; }
1295 0 : bool supportsElemSideQpArg() const override final { return true; }
1296 :
1297 : private:
1298 17980204 : ValueType evaluate(const ElemArg &, const StateArg &) const override { return _value; }
1299 2285783 : ValueType evaluate(const FaceArg &, const StateArg &) const override { return _value; }
1300 3977747 : ValueType evaluate(const ElemQpArg &, const StateArg &) const override { return _value; }
1301 270366 : ValueType evaluate(const ElemSideQpArg &, const StateArg &) const override { return _value; }
1302 196 : ValueType evaluate(const ElemPointArg &, const StateArg &) const override { return _value; }
1303 45415 : ValueType evaluate(const NodeArg &, const StateArg &) const override { return _value; }
1304 :
1305 10 : GradientType evaluateGradient(const ElemArg &, const StateArg &) const override { return 0; }
1306 10 : GradientType evaluateGradient(const FaceArg &, const StateArg &) const override { return 0; }
1307 10 : GradientType evaluateGradient(const ElemQpArg &, const StateArg &) const override { return 0; }
1308 10 : GradientType evaluateGradient(const ElemSideQpArg &, const StateArg &) const override
1309 : {
1310 10 : return 0;
1311 : }
1312 10 : GradientType evaluateGradient(const ElemPointArg &, const StateArg &) const override { return 0; }
1313 10 : GradientType evaluateGradient(const NodeArg &, const StateArg &) const override { return 0; }
1314 :
1315 10 : DotType evaluateDot(const ElemArg &, const StateArg &) const override { return 0; }
1316 10 : DotType evaluateDot(const FaceArg &, const StateArg &) const override { return 0; }
1317 10 : DotType evaluateDot(const ElemQpArg &, const StateArg &) const override { return 0; }
1318 10 : DotType evaluateDot(const ElemSideQpArg &, const StateArg &) const override { return 0; }
1319 10 : DotType evaluateDot(const ElemPointArg &, const StateArg &) const override { return 0; }
1320 10 : DotType evaluateDot(const NodeArg &, const StateArg &) const override { return 0; }
1321 :
1322 10 : GradientType evaluateGradDot(const ElemArg &, const StateArg &) const override { return 0; }
1323 10 : GradientType evaluateGradDot(const FaceArg &, const StateArg &) const override { return 0; }
1324 10 : GradientType evaluateGradDot(const ElemQpArg &, const StateArg &) const override { return 0; }
1325 10 : GradientType evaluateGradDot(const ElemSideQpArg &, const StateArg &) const override { return 0; }
1326 10 : GradientType evaluateGradDot(const ElemPointArg &, const StateArg &) const override { return 0; }
1327 10 : GradientType evaluateGradDot(const NodeArg &, const StateArg &) const override { return 0; }
1328 :
1329 : private:
1330 : ValueType _value;
1331 : };
1332 :
1333 : /**
1334 : * A functor that serves as a placeholder during the simulation setup phase if a functor consumer
1335 : * requests a functor that has not yet been constructed.
1336 : */
1337 : template <typename T>
1338 : class NullFunctor final : public FunctorBase<T>
1339 : {
1340 : public:
1341 : using typename FunctorBase<T>::FunctorType;
1342 : using typename FunctorBase<T>::ValueType;
1343 : using typename FunctorBase<T>::GradientType;
1344 : using typename FunctorBase<T>::DotType;
1345 :
1346 1088 : NullFunctor() : FunctorBase<T>("null") {}
1347 :
1348 : // For backwards compatiblity of unit testing
1349 : bool hasFaceSide(const FaceInfo & fi, bool) const override;
1350 :
1351 0 : bool supportsFaceArg() const override final { return false; }
1352 0 : bool supportsElemSideQpArg() const override final { return false; }
1353 :
1354 : private:
1355 2 : ValueType evaluate(const ElemArg &, const StateArg &) const override
1356 : {
1357 2 : mooseError("We should never get here. If you have, contact a MOOSE developer and tell them "
1358 : "they've written broken code");
1359 : }
1360 2 : ValueType evaluate(const FaceArg &, const StateArg &) const override
1361 : {
1362 2 : mooseError("We should never get here. If you have, contact a MOOSE developer and tell them "
1363 : "they've written broken code");
1364 : }
1365 2 : ValueType evaluate(const ElemQpArg &, const StateArg &) const override
1366 : {
1367 2 : mooseError("We should never get here. If you have, contact a MOOSE developer and tell them "
1368 : "they've written broken code");
1369 : }
1370 2 : ValueType evaluate(const ElemSideQpArg &, const StateArg &) const override
1371 : {
1372 2 : mooseError("We should never get here. If you have, contact a MOOSE developer and tell them "
1373 : "they've written broken code");
1374 : }
1375 2 : ValueType evaluate(const ElemPointArg &, const StateArg &) const override
1376 : {
1377 2 : mooseError("We should never get here. If you have, contact a MOOSE developer and tell them "
1378 : "they've written broken code");
1379 : }
1380 2 : ValueType evaluate(const NodeArg &, const StateArg &) const override
1381 : {
1382 2 : mooseError("We should never get here. If you have, contact a MOOSE developer and tell them "
1383 : "they've written broken code");
1384 : }
1385 : };
1386 :
1387 : template <typename T>
1388 : bool
1389 0 : NullFunctor<T>::hasFaceSide(const FaceInfo &, const bool) const
1390 : {
1391 : // For backwards compatiblity of unit testing
1392 0 : return true;
1393 : }
1394 : }
|