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 "MooseVariableBase.h" 13 : #include "MooseFunctor.h" 14 : #include "SystemBase.h" 15 : 16 : // libMesh forward declarations 17 : namespace libMesh 18 : { 19 : template <typename T> 20 : class NumericVector; 21 : } 22 : 23 : class Assembly; 24 : class TimeIntegrator; 25 : 26 : /** 27 : * Class for scalar variables (they are different). 28 : */ 29 : class MooseVariableScalar : public MooseVariableBase, public Moose::FunctorBase<ADReal> 30 : { 31 : public: 32 : static InputParameters validParams(); 33 : 34 : MooseVariableScalar(const InputParameters & parameters); 35 : virtual ~MooseVariableScalar(); 36 : 37 : /** 38 : * Fill out the VariableValue arrays from the system solution vector 39 : * @param reinit_for_derivative_reordering A flag indicating whether we are reinitializing for the 40 : * purpose of re-ordering derivative information for ADNodalBCs 41 : */ 42 : void reinit(bool reinit_for_derivative_reordering = false); 43 : 44 66716 : const VariableValue & sln() const { return _u; } 45 : 46 : /** 47 : * Return the solution with derivative information 48 : */ 49 : const ADVariableValue & adSln() const; 50 : 51 : const VariableValue & slnOld() const; 52 : const VariableValue & slnOlder() const; 53 : const VariableValue & vectorTagSln(TagID tag) const; 54 : const VariableValue & matrixTagSln(TagID tag) const; 55 : 56 : const VariableValue & uDot() const; 57 : const VariableValue & uDotDot() const; 58 : const VariableValue & uDotOld() const; 59 : const VariableValue & uDotDotOld() const; 60 : const VariableValue & duDotDu() const; 61 : const VariableValue & duDotDotDu() const; 62 : 63 : /** 64 : * Return the first derivative of the solution with derivative information 65 : */ 66 : const ADVariableValue & adUDot() const; 67 : 68 : /** 69 : * Set the nodal value for this variable (to keep everything up to date 70 : */ 71 : void setValue(unsigned int i, Number value); 72 : 73 : /** 74 : * Set all of the values of this scalar variable to the same value 75 : */ 76 : void setValues(Number value); 77 : 78 : void insert(NumericVector<Number> & soln); 79 : 80 : /** 81 : * The oldest solution state that is requested for this variable 82 : * (0 = current, 1 = old, 2 = older, etc). 83 : */ 84 : unsigned int oldestSolutionStateRequested() const; 85 : 86 123006 : void setActiveTags(const std::set<TagID> & vtags) override { _required_vector_tags = vtags; } 87 : 88 : virtual void sizeMatrixTagData() override; 89 : 90 0 : bool supportsFaceArg() const override final { return true; } 91 0 : bool supportsElemSideQpArg() const override final { return true; } 92 : 93 : protected: 94 : /// The value of scalar variable 95 : VariableValue _u; 96 : /// The old value of scalar variable 97 : VariableValue _u_old; 98 : /// The older value of scalar variable 99 : VariableValue _u_older; 100 : /// Tagged vectors 101 : std::vector<VariableValue> _vector_tag_u; 102 : /// Only cache data when need it 103 : mutable std::vector<bool> _need_vector_tag_u; 104 : /// Tagged matrices 105 : std::vector<VariableValue> _matrix_tag_u; 106 : /// Only cache data when need it 107 : mutable std::vector<bool> _need_matrix_tag_u; 108 : 109 : VariableValue _u_dot; 110 : VariableValue _u_dotdot; 111 : VariableValue _u_dot_old; 112 : VariableValue _u_dotdot_old; 113 : VariableValue _du_dot_du; 114 : VariableValue _du_dotdot_du; 115 : 116 : mutable bool _need_u_dot; 117 : mutable bool _need_u_dotdot; 118 : mutable bool _need_u_dot_old; 119 : mutable bool _need_u_dotdot_old; 120 : mutable bool _need_du_dot_du; 121 : mutable bool _need_du_dotdot_du; 122 : 123 : /// Whether or not the old solution is needed 124 : mutable bool _need_u_old; 125 : /// Whether or not the older solution is needed 126 : mutable bool _need_u_older; 127 : 128 : /// Whether any AD calculations are needed 129 : mutable bool _need_ad; 130 : /// whether ad_u is needed 131 : mutable bool _need_ad_u; 132 : /// whether ad_u_dot is needed 133 : mutable bool _need_ad_u_dot; 134 : /// The scalar solution with derivative information 135 : ADVariableValue _ad_u; 136 : /// The first derivative of the scalar solution with derivative information 137 : ADVariableValue _ad_u_dot; 138 : 139 : private: 140 : using typename Moose::FunctorBase<ADReal>::ValueType; 141 : using typename Moose::FunctorBase<ADReal>::GradientType; 142 : using typename Moose::FunctorBase<ADReal>::DotType; 143 : 144 : using ElemArg = Moose::ElemArg; 145 : using ElemQpArg = Moose::ElemQpArg; 146 : using ElemSideQpArg = Moose::ElemSideQpArg; 147 : using FaceArg = Moose::FaceArg; 148 : using ElemPointArg = Moose::ElemPointArg; 149 : using NodeArg = Moose::NodeArg; 150 : 151 : ValueType evaluate(const ElemArg & elem, const Moose::StateArg & state) const override final; 152 : ValueType evaluate(const FaceArg & face, const Moose::StateArg & state) const override final; 153 : ValueType evaluate(const ElemQpArg & qp, const Moose::StateArg & state) const override final; 154 : ValueType evaluate(const ElemSideQpArg & elem_side_qp, 155 : const Moose::StateArg & state) const override final; 156 : ValueType evaluate(const ElemPointArg & elem_point, 157 : const Moose::StateArg & state) const override final; 158 : ValueType evaluate(const NodeArg & node, const Moose::StateArg & state) const override final; 159 : ValueType evaluate(const Moose::StateArg & state) const; 160 : 161 : GradientType evaluateGradient(const ElemArg & elem, 162 : const Moose::StateArg & state) const override final; 163 : GradientType evaluateGradient(const FaceArg & face, 164 : const Moose::StateArg & state) const override final; 165 : GradientType evaluateGradient(const ElemQpArg & qp, 166 : const Moose::StateArg & state) const override final; 167 : GradientType evaluateGradient(const ElemSideQpArg & elem_side_qp, 168 : const Moose::StateArg & state) const override final; 169 : GradientType evaluateGradient(const ElemPointArg & elem_point, 170 : const Moose::StateArg & state) const override final; 171 : GradientType evaluateGradient(const NodeArg & node, 172 : const Moose::StateArg & state) const override final; 173 : 174 : /** 175 : * Adds derivative information to the scalar variable value arrays 176 : * @param nodal_ordering Whether we are doing a nodal ordering of the derivative vector, e.g. 177 : * whether the spacing between variable groups in the derivative vector is equal to the 178 : * number of dofs per node or the number of dofs per elem 179 : */ 180 : void computeAD(bool nodal_ordering); 181 : 182 : /// The set of vector tags we need to evaluate 183 : std::set<TagID> _required_vector_tags; 184 : };