https://mooseframework.inl.gov
RawValueFunctor.h
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 #pragma once
11 
12 #include "MooseFunctor.h"
13 
14 namespace Moose
15 {
16 template <typename T>
17 class RawValueFunctor : public FunctorBase<T>
18 {
19 public:
20  using typename Moose::FunctorBase<T>::ValueType;
22  using typename Moose::FunctorBase<T>::DotType;
23 
24  RawValueFunctor(const FunctorBase<typename ADType<T>::type> & ad_functor)
25  : FunctorBase<T>(ad_functor.functorName() + "_raw_value", {EXEC_ALWAYS}),
26  _ad_functor(ad_functor)
27  {
28  }
29 
30  virtual bool isExtrapolatedBoundaryFace(const FaceInfo & fi,
31  const Elem * const elem,
32  const Moose::StateArg & state) const override
33  {
34  return _ad_functor.isExtrapolatedBoundaryFace(fi, elem, state);
35  }
36  virtual bool isConstant() const override { return _ad_functor.isConstant(); }
37  virtual bool hasBlocks(const SubdomainID id) const override { return _ad_functor.hasBlocks(id); }
38  virtual bool hasFaceSide(const FaceInfo & fi, const bool fi_elem_side) const override
39  {
40  return _ad_functor.hasFaceSide(fi, fi_elem_side);
41  }
42 
43  bool supportsFaceArg() const override final { return _ad_functor.supportsFaceArg(); }
44  bool supportsElemSideQpArg() const override final { return _ad_functor.supportsElemSideQpArg(); }
45 
46 protected:
48 
51  ValueType evaluate(const ElemArg & elem, const StateArg & state) const override
52  {
53  return MetaPhysicL::raw_value(_ad_functor(elem, state));
54  }
55  ValueType evaluate(const FaceArg & face, const StateArg & state) const override
56  {
57  return MetaPhysicL::raw_value(_ad_functor(face, state));
58  }
59  ValueType evaluate(const ElemQpArg & qp, const StateArg & state) const override
60  {
61  return MetaPhysicL::raw_value(_ad_functor(qp, state));
62  }
63  ValueType evaluate(const ElemSideQpArg & qp, const StateArg & state) const override
64  {
65  return MetaPhysicL::raw_value(_ad_functor(qp, state));
66  }
67  ValueType evaluate(const ElemPointArg & elem_point, const StateArg & state) const override
68  {
69  return MetaPhysicL::raw_value(_ad_functor(elem_point, state));
70  }
71  ValueType evaluate(const NodeArg & node, const StateArg & state) const override
72  {
73  return MetaPhysicL::raw_value(_ad_functor(node, state));
74  }
75 
76  GradientType evaluateGradient(const ElemArg & elem, const StateArg & state) const override
77  {
78  return MetaPhysicL::raw_value(_ad_functor.gradient(elem, state));
79  }
80  GradientType evaluateGradient(const FaceArg & face, const StateArg & state) const override
81  {
82  return MetaPhysicL::raw_value(_ad_functor.gradient(face, state));
83  }
84  GradientType evaluateGradient(const ElemQpArg & qp, const StateArg & state) const override
85  {
86  return MetaPhysicL::raw_value(_ad_functor.gradient(qp, state));
87  }
88  GradientType evaluateGradient(const ElemSideQpArg & qp, const StateArg & state) const override
89  {
90  return MetaPhysicL::raw_value(_ad_functor.gradient(qp, state));
91  }
93  const StateArg & state) const override
94  {
95  return MetaPhysicL::raw_value(_ad_functor.gradient(elem_point, state));
96  }
97  GradientType evaluateGradient(const NodeArg & node, const StateArg & state) const override
98  {
99  return MetaPhysicL::raw_value(_ad_functor.gradient(node, state));
100  }
101 
102  DotType evaluateDot(const ElemArg & elem, const StateArg & state) const override
103  {
104  return MetaPhysicL::raw_value(_ad_functor.dot(elem, state));
105  }
106  DotType evaluateDot(const FaceArg & face, const StateArg & state) const override
107  {
108  return MetaPhysicL::raw_value(_ad_functor.dot(face, state));
109  }
110  DotType evaluateDot(const ElemQpArg & qp, const StateArg & state) const override
111  {
112  return MetaPhysicL::raw_value(_ad_functor.dot(qp, state));
113  }
114  DotType evaluateDot(const ElemSideQpArg & qp, const StateArg & state) const override
115  {
116  return MetaPhysicL::raw_value(_ad_functor.dot(qp, state));
117  }
118  DotType evaluateDot(const ElemPointArg & elem_point, const StateArg & state) const override
119  {
120  return MetaPhysicL::raw_value(_ad_functor.dot(elem_point, state));
121  }
122  DotType evaluateDot(const NodeArg & node, const StateArg & state) const override
123  {
124  return MetaPhysicL::raw_value(_ad_functor.dot(node, state));
125  }
126 
127  GradientType evaluateGradDot(const ElemArg & elem, const StateArg & state) const override
128  {
129  return MetaPhysicL::raw_value(_ad_functor.gradDot(elem, state));
130  }
131  GradientType evaluateGradDot(const FaceArg & face, const StateArg & state) const override
132  {
133  return MetaPhysicL::raw_value(_ad_functor.gradDot(face, state));
134  }
135  GradientType evaluateGradDot(const ElemQpArg & qp, const StateArg & state) const override
136  {
137  return MetaPhysicL::raw_value(_ad_functor.gradDot(qp, state));
138  }
139  GradientType evaluateGradDot(const ElemSideQpArg & qp, const StateArg & state) const override
140  {
141  return MetaPhysicL::raw_value(_ad_functor.gradDot(qp, state));
142  }
144  const StateArg & state) const override
145  {
146  return MetaPhysicL::raw_value(_ad_functor.gradDot(elem_point, state));
147  }
148  GradientType evaluateGradDot(const NodeArg & node, const StateArg & state) const override
149  {
150  return MetaPhysicL::raw_value(_ad_functor.gradDot(node, state));
151  }
153 
154 private:
157 };
158 }
virtual bool isExtrapolatedBoundaryFace(const FaceInfo &fi, const Elem *const elem, const Moose::StateArg &state) const override
Returns whether this (sided) face is an extrapolated boundary face for this functor.
const FunctorBase< typename ADType< T >::type > & _ad_functor
Our wrapped AD object.
DotType evaluateDot(const ElemSideQpArg &qp, const StateArg &state) const override
DotType evaluateDot(const ElemPointArg &elem_point, const StateArg &state) const override
Evaluate the functor time derivative with a given element and point.
Base class template for functor objects.
Definition: MooseFunctor.h:137
DotType evaluateDot(const NodeArg &node, const StateArg &state) const override
ValueType evaluate(const ElemQpArg &qp, const StateArg &state) const override
GradientType evaluateGradDot(const NodeArg &node, const StateArg &state) const override
virtual bool hasFaceSide(const FaceInfo &fi, const bool fi_elem_side) const override
GradientType evaluateGradDot(const ElemPointArg &elem_point, const StateArg &state) const override
Evaluate the functor gradient-dot with a given element and point.
GradientType evaluateGradient(const NodeArg &node, const StateArg &state) const override
GradientType evaluateGradDot(const ElemArg &elem, const StateArg &state) const override
Evaluate the functor gradient-dot with a given element.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
GradientType evaluateGradient(const ElemArg &elem, const StateArg &state) const override
Evaluate the functor gradient with a given element.
virtual bool hasBlocks(const SubdomainID id) const override
Returns whether the functor is defined on this block.
A structure that is used to evaluate Moose functors at an arbitrary physical point contained within a...
GradientType evaluateGradient(const ElemSideQpArg &qp, const StateArg &state) const override
const ExecFlagType EXEC_ALWAYS
Definition: Moose.C:47
ValueType evaluate(const ElemPointArg &elem_point, const StateArg &state) const override
Evaluate the functor with a given element and point.
typename FunctorReturnType< T, FunctorEvaluationKind::Gradient >::type GradientType
This rigmarole makes it so that a user can create functors that return containers (std::vector...
Definition: MooseFunctor.h:149
virtual bool isConstant() const override
Returns true if this functor is a constant.
GradientType evaluateGradDot(const ElemSideQpArg &qp, const StateArg &state) const override
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
GradientType evaluateGradient(const FaceArg &face, const StateArg &state) const override
A structure defining a "face" evaluation calling argument for Moose functors.
ValueType evaluate(const NodeArg &node, const StateArg &state) const override
GradientType evaluateGradDot(const ElemQpArg &qp, const StateArg &state) const override
GradientType evaluateGradient(const ElemPointArg &elem_point, const StateArg &state) const override
Evaluate the functor gradient with a given element and point.
RawValueFunctor(const FunctorBase< typename ADType< T >::type > &ad_functor)
A structure that is used to evaluate Moose functors logically at an element/cell center.
ValueType evaluate(const ElemArg &elem, const StateArg &state) const override
Forward calls to wrapped object.
Argument for requesting functor evaluation at a quadrature point location in an element.
bool supportsElemSideQpArg() const override final
Whether this functor supports evaluation with ElemSideQpArg.
const MooseFunctorName & functorName() const
Return the functor name.
Definition: MooseFunctor.h:170
ValueType evaluate(const ElemSideQpArg &qp, const StateArg &state) const override
GradientType evaluateGradient(const ElemQpArg &qp, const StateArg &state) const override
DotType evaluateDot(const ElemQpArg &qp, const StateArg &state) const override
State argument for evaluating functors.
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
DotType evaluateDot(const FaceArg &face, const StateArg &state) const override
GradientType evaluateGradDot(const FaceArg &face, const StateArg &state) const override
DotType evaluateDot(const ElemArg &elem, const StateArg &state) const override
Evaluate the functor time derivative with a given element.
Argument for requesting functor evaluation at quadrature point locations on an element side...
ValueType evaluate(const FaceArg &face, const StateArg &state) const override
bool supportsFaceArg() const override final
Whether this functor supports evaluation with FaceArg.