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 "MooseTypes.h" 13 : #include "MooseVariableField.h" 14 : #include "SubProblem.h" 15 : #include "MooseMesh.h" 16 : #include "MooseVariableDataLinearFV.h" 17 : #include "GradientLimiterType.h" 18 : 19 : #include "libmesh/numeric_vector.h" 20 : #include "libmesh/dof_map.h" 21 : #include "libmesh/elem.h" 22 : #include "libmesh/quadrature.h" 23 : #include "libmesh/dense_vector.h" 24 : #include "libmesh/enum_fe_family.h" 25 : 26 : template <typename> 27 : class MooseLinearVariableFV; 28 : 29 : typedef MooseLinearVariableFV<Real> MooseLinearVariableFVReal; 30 : class AuxiliarySystem; 31 : class FVDirichletBCBase; 32 : class FVFluxBC; 33 : class LinearFVBoundaryCondition; 34 : class LinearSystem; 35 : 36 : namespace libMesh 37 : { 38 : template <typename> 39 : class NumericVector; 40 : } 41 : 42 : /// This class provides variable solution interface for linear 43 : /// finite volume problems. 44 : /// This class is designed to store gradient information when enabled. 45 : template <typename OutputType> 46 : class MooseLinearVariableFV : public MooseVariableField<OutputType> 47 : { 48 : public: 49 : using OutputGradient = typename MooseVariableField<OutputType>::OutputGradient; 50 : using OutputSecond = typename MooseVariableField<OutputType>::OutputSecond; 51 : using OutputDivergence = typename MooseVariableField<OutputType>::OutputDivergence; 52 : 53 : using FieldVariableValue = typename MooseVariableField<OutputType>::FieldVariableValue; 54 : using FieldVariableGradient = typename MooseVariableField<OutputType>::FieldVariableGradient; 55 : using FieldVariableSecond = typename MooseVariableField<OutputType>::FieldVariableSecond; 56 : using FieldVariableCurl = typename MooseVariableField<OutputType>::FieldVariableCurl; 57 : using FieldVariableDivergence = typename MooseVariableField<OutputType>::FieldVariableDivergence; 58 : 59 : using OutputShape = typename MooseVariableField<OutputType>::OutputShape; 60 : using OutputShapeGradient = typename MooseVariableField<OutputType>::OutputShapeGradient; 61 : using OutputShapeSecond = typename MooseVariableField<OutputType>::OutputShapeSecond; 62 : using OutputShapeDivergence = typename MooseVariableField<OutputType>::OutputShapeDivergence; 63 : 64 : using DofValue = typename MooseVariableField<OutputType>::DofValue; 65 : using DofValues = typename MooseVariableField<OutputType>::DofValues; 66 : using ADDofValue = typename MooseVariableField<OutputType>::ADDofValue; 67 : using ADDofValues = typename MooseVariableField<OutputType>::ADDofValues; 68 : 69 : using FieldVariablePhiValue = typename MooseVariableField<OutputType>::FieldVariablePhiValue; 70 : using FieldVariablePhiGradient = 71 : typename MooseVariableField<OutputType>::FieldVariablePhiGradient; 72 : using FieldVariablePhiSecond = typename MooseVariableField<OutputType>::FieldVariablePhiSecond; 73 : using FieldVariablePhiDivergence = 74 : typename MooseVariableField<OutputType>::FieldVariablePhiDivergence; 75 : using ElemQpArg = Moose::ElemQpArg; 76 : using ElemSideQpArg = Moose::ElemSideQpArg; 77 : using ElemArg = Moose::ElemArg; 78 : using FaceArg = Moose::FaceArg; 79 : using StateArg = Moose::StateArg; 80 : using NodeArg = Moose::NodeArg; 81 : using ElemPointArg = Moose::ElemPointArg; 82 : using typename MooseVariableField<OutputType>::ValueType; 83 : using typename MooseVariableField<OutputType>::DotType; 84 : using typename MooseVariableField<OutputType>::GradientType; 85 : 86 : static InputParameters validParams(); 87 : MooseLinearVariableFV(const InputParameters & parameters); 88 : 89 1101776 : virtual bool isFV() const override { return true; } 90 : 91 : /** 92 : * If the variable has a dirichlet boundary condition at face described by \p fi . 93 : */ 94 : virtual bool isDirichletBoundaryFace(const FaceInfo & fi) const; 95 : 96 : /** 97 : * Switch to request cell gradient computations. 98 : */ 99 1190 : void computeCellGradients() { _needs_cell_gradients = true; } 100 : 101 : /** 102 : * Switch to request cell gradient computations with an optional gradient limiter. 103 : * 104 : * `GradientLimiterType::None` is equivalent to requesting the regular gradients only. 105 : */ 106 : void computeCellGradients(const Moose::FV::GradientLimiterType limiter_type); 107 : 108 : /** 109 : * Switch to request limited cell gradient computations. 110 : * 111 : * Limited gradients are stored in limiter-specific containers on the system and are computed 112 : * using the raw cell gradients. 113 : */ 114 : void computeCellLimitedGradients(const Moose::FV::GradientLimiterType limiter_type); 115 : 116 : /** 117 : * Check if cell gradient computations were requested for this variable. 118 : */ 119 100323 : virtual bool needsGradientVectorStorage() const override { return _needs_cell_gradients; } 120 : 121 : virtual bool isExtrapolatedBoundaryFace(const FaceInfo & fi, 122 : const Elem * elem, 123 : const Moose::StateArg & state) const override; 124 : 125 : /** 126 : * Get the variable gradient at a cell center. 127 : * @param elem_info The ElemInfo of the cell where we need the gradient 128 : * @param state State argument describing which solution state to evaluate 129 : */ 130 : VectorValue<Real> gradSln(const ElemInfo & elem_info, const StateArg & state) const; 131 : 132 : /** 133 : * Get one raw gradient component at a cell center without materializing the full gradient. 134 : * @param elem_info The ElemInfo of the cell where we need the gradient 135 : * @param component The gradient component to retrieve 136 : */ 137 : Real gradSlnComponent(const ElemInfo & elem_info, unsigned int component) const; 138 : 139 : /** 140 : * Get either the raw or limited gradient at a cell center. 141 : * @param elem_info The ElemInfo of the cell where we need the gradient 142 : * @param state State argument describing which solution state to evaluate 143 : * @param limiter_type The limiter type used to compute/store limited gradients 144 : */ 145 : VectorValue<Real> gradSln(const ElemInfo & elem_info, 146 : const StateArg & state, 147 : const Moose::FV::GradientLimiterType limiter_type) const; 148 : 149 : /** 150 : * Get the limited gradient at a cell center. 151 : * @param elem_info The ElemInfo of the cell where we need the gradient 152 : * @param state State argument describing which solution state to evaluate 153 : * @param limiter_type The limiter type used to compute/store limited gradients 154 : */ 155 : VectorValue<Real> limitedGradSln(const ElemInfo & elem_info, 156 : const StateArg & state, 157 : const Moose::FV::GradientLimiterType limiter_type) const; 158 : 159 : /** 160 : * Compute interpolated gradient on the provided face. 161 : * @param fi The face for which to retrieve the gradient 162 : * @param state State argument describing which solution state to evaluate 163 : */ 164 : VectorValue<Real> gradSln(const FaceInfo & fi, const StateArg & state) const; 165 : 166 : /** 167 : * Compute interpolated raw/limited gradient on the provided face. 168 : * @param fi The face for which to retrieve the gradient 169 : * @param state State argument describing which solution state to evaluate 170 : * @param limiter_type The limiter type used to compute/store limited gradients 171 : */ 172 : VectorValue<Real> gradSln(const FaceInfo & fi, 173 : const StateArg & state, 174 : const Moose::FV::GradientLimiterType limiter_type) const; 175 : 176 : /** 177 : * Compute interpolated limited gradient on the provided face. 178 : * @param fi The face for which to retrieve the gradient 179 : * @param state State argument describing which solution state to evaluate 180 : * @param limiter_type The limiter type used to compute/store limited gradients 181 : */ 182 : VectorValue<Real> limitedGradSln(const FaceInfo & fi, 183 : const StateArg & state, 184 : const Moose::FV::GradientLimiterType limiter_type) const; 185 : 186 : virtual void initialSetup() override; 187 : virtual void timestepSetup() override; 188 : 189 : /** 190 : * Get the solution value for the provided element and seed the derivative for the corresponding 191 : * dof index 192 : * @param elem_info The element to retrieve the solution value for 193 : * @param state State argument which describes at what time / solution iteration state we want to 194 : * evaluate the variable 195 : */ 196 : Real getElemValue(const ElemInfo & elem_info, const StateArg & state) const; 197 : 198 : /** 199 : * Get the boundary condition object which corresponds to the given boundary ID 200 : * @param bd_id The boundary ID whose condition should be fetched 201 : */ 202 : LinearFVBoundaryCondition * getBoundaryCondition(const BoundaryID bd_id) const; 203 : 204 1164 : const std::unordered_map<BoundaryID, LinearFVBoundaryCondition *> & getBoundaryConditionMap() 205 : { 206 1164 : return _boundary_id_to_bc; 207 : } 208 : 209 396310 : virtual void prepareIC() override {} 210 : 211 1902 : virtual bool isNodal() const override final { return false; } 212 : 213 0 : virtual bool hasDoFsOnNodes() const override final { return false; } 214 : 215 0 : virtual bool isNodalDefined() const override final { return false; } 216 : 217 0 : virtual bool supportsFaceArg() const override final { return true; } 218 0 : virtual bool supportsElemSideQpArg() const override final { return false; } 219 : 220 : virtual const Elem * const & currentElem() const override; 221 : 222 0 : virtual bool computingSecond() const override final { return false; } 223 0 : virtual bool computingCurl() const override final { return false; } 224 0 : virtual bool computingDiv() const override final { return false; } 225 0 : virtual bool usesSecondPhiNeighbor() const override final { return false; } 226 : 227 : virtual void sizeMatrixTagData() override; 228 : 229 : protected: 230 : /// Throw an error when somebody requests time-related data from this variable 231 : [[noreturn]] void timeIntegratorError() const; 232 : 233 : /// Throw and error when somebody requests lower-dimensional data from this variable 234 : [[noreturn]] void lowerDError() const; 235 : 236 : /// Throw an error when somebody wants to use this variable as a nodal variable 237 : [[noreturn]] void nodalError() const; 238 : 239 : /// Throw an error when somebody wants to use this variable with automatic differentiation 240 : [[noreturn]] void adError() const; 241 : 242 : /// Throw an error when somebody requests gradients at a non-current solution state 243 : [[noreturn]] void gradientStateError(const StateArg & state) const; 244 : 245 : /** 246 : * Setup the boundary to Dirichlet BC map 247 : */ 248 : void cacheBoundaryBCMap(); 249 : 250 : usingMooseVariableBaseMembers; 251 : 252 : /// Boolean to check if this variable needs gradient computations. 253 : bool _needs_cell_gradients; 254 : 255 : /// Temporary storage for the cell gradient to avoid unnecessary allocations. 256 : mutable RealVectorValue _cell_gradient; 257 : 258 : /// Owning concrete system pointers. One will be null. 259 : LinearSystem * const _linear_system; 260 : AuxiliarySystem * const _auxiliary_system; 261 : 262 : /// Pointer to the unlimited cell gradient stored by the owning concrete system 263 : const std::vector<std::unique_ptr<libMesh::NumericVector<libMesh::Number>>> & _grad_container; 264 : 265 : /// Holder for all the data associated with the "main" element. The data in this is 266 : /// mainly used by finite element-based loops such as the postprocessor and auxkernel 267 : /// loops 268 : std::unique_ptr<MooseVariableDataLinearFV<OutputType>> _element_data; 269 : 270 : /// Holder for all the data associated with the "neighbor" element. The data in this is 271 : /// mainly used by finite element-based loops such as the postprocessor and auxkernel 272 : /// loops 273 : std::unique_ptr<MooseVariableDataLinearFV<OutputType>> _neighbor_data; 274 : 275 : /// Map for easily accessing the boundary conditions based on the boundary IDs. 276 : /// We assume that each boundary has one boundary condition only. 277 : std::unordered_map<BoundaryID, LinearFVBoundaryCondition *> _boundary_id_to_bc; 278 : 279 : /// Cache the number of the system this variable belongs to 280 : const unsigned int _sys_num; 281 : 282 : friend void Moose::initDofIndices<>(MooseLinearVariableFV<OutputType> &, const Elem &); 283 : 284 : private: 285 : using MooseVariableField<OutputType>::evaluate; 286 : using MooseVariableField<OutputType>::evaluateGradient; 287 : using MooseVariableField<OutputType>::evaluateDot; 288 : 289 : virtual ValueType evaluate(const ElemArg & elem, const StateArg &) const override final; 290 : virtual ValueType evaluate(const FaceArg & face, const StateArg &) const override final; 291 : virtual ValueType evaluate(const NodeArg & node, const StateArg &) const override final; 292 : virtual ValueType evaluate(const ElemPointArg & elem_point, 293 : const StateArg & state) const override final; 294 : virtual ValueType evaluate(const ElemQpArg & elem_qp, 295 : const StateArg & state) const override final; 296 : virtual ValueType evaluate(const ElemSideQpArg & elem_side_qp, 297 : const StateArg & state) const override final; 298 : virtual GradientType evaluateGradient(const ElemQpArg & qp_arg, 299 : const StateArg &) const override final; 300 : virtual GradientType evaluateGradient(const ElemArg & elem_arg, 301 : const StateArg &) const override final; 302 : virtual GradientType evaluateGradient(const FaceArg & face, 303 : const StateArg &) const override final; 304 : virtual DotType evaluateDot(const ElemArg & elem, const StateArg &) const override final; 305 : 306 : /// The current (ghosted) solution. Note that this needs to be stored as a reference to a pointer 307 : /// because the solution might not exist at the time that this variable is constructed, so we 308 : /// cannot safely dereference at that time 309 : const libMesh::NumericVector<libMesh::Number> * const & _solution; 310 : 311 : /// Shape functions, only used when we are postprocessing or using this variable 312 : /// in an auxiliary system 313 : const FieldVariablePhiValue & _phi; 314 : const FieldVariablePhiGradient & _grad_phi; 315 : const FieldVariablePhiValue & _phi_face; 316 : const FieldVariablePhiGradient & _grad_phi_face; 317 : const FieldVariablePhiValue & _phi_face_neighbor; 318 : const FieldVariablePhiGradient & _grad_phi_face_neighbor; 319 : const FieldVariablePhiValue & _phi_neighbor; 320 : const FieldVariablePhiGradient & _grad_phi_neighbor; 321 : 322 : public: 323 : // ********************************************************************************* 324 : // ********************************************************************************* 325 : // The following functions are separated here because they are not essential for the 326 : // solver but are necessary to interface with the auxiliary and postprocessor 327 : // systems. 328 : // ********************************************************************************* 329 : // ********************************************************************************* 330 : 331 : virtual void setDofValue(const DofValue & /*value*/, unsigned int /*index*/) override; 332 : 333 : virtual void getDofIndices(const Elem * elem, 334 : std::vector<dof_id_type> & dof_indices) const override; 335 : 336 : virtual void setDofValues(const DenseVector<DofValue> & values) override; 337 : 338 : virtual void clearDofIndices() override; 339 : 340 434 : virtual unsigned int numberOfDofs() const override final { return 1; } 341 0 : virtual unsigned int numberOfDofsNeighbor() override final { return 1; } 342 : 343 : virtual unsigned int oldestSolutionStateRequested() const override final; 344 : 345 : virtual void clearAllDofIndices() override final; 346 : 347 : [[noreturn]] virtual const std::vector<dof_id_type> & dofIndicesLower() const override final; 348 : [[noreturn]] virtual const FieldVariablePhiValue & phiLower() const override; 349 : 350 : // Overriding these to make sure nothing happens during residual/jacobian setup. 351 : // The only time this can actually happen is when residual setup is called on the auxiliary 352 : // system. 353 20 : virtual void residualSetup() override {} 354 290 : virtual void jacobianSetup() override {} 355 : 356 0 : virtual libMesh::FEContinuity getContinuity() const override 357 : { 358 0 : return _element_data->getContinuity(); 359 : }; 360 : 361 : virtual void setNodalValue(const OutputType & value, unsigned int idx = 0) override; 362 : 363 : [[noreturn]] virtual const DofValues & nodalVectorTagValue(TagID) const override; 364 : 365 : virtual const std::vector<dof_id_type> & dofIndices() const final; 366 : virtual const std::vector<dof_id_type> & dofIndicesNeighbor() const final; 367 : 368 575365 : virtual void prepare() override final {} 369 0 : virtual void prepareNeighbor() override final {} 370 0 : virtual void prepareAux() override final {} 371 0 : virtual void reinitNode() override final {} 372 0 : virtual void reinitNodes(const std::vector<dof_id_type> & /*nodes*/) override final {} 373 0 : virtual void reinitNodesNeighbor(const std::vector<dof_id_type> & /*nodes*/) override final {} 374 1549 : virtual void reinitAux() override final {} 375 0 : virtual void reinitAuxNeighbor() override final {} 376 0 : virtual void prepareLowerD() override final {} 377 : 378 : virtual void computeElemValuesFace() override; 379 : virtual void computeNeighborValuesFace() override; 380 : virtual void computeNeighborValues() override; 381 : virtual void computeLowerDValues() override final; 382 : 383 : virtual void computeNodalNeighborValues() override final; 384 : virtual void computeNodalValues() override final; 385 : 386 : virtual void computeElemValues() override; 387 0 : virtual void computeFaceValues(const FaceInfo & /*fi*/) override {} 388 : 389 : virtual void setLowerDofValues(const DenseVector<DofValue> & values) override; 390 : 391 : virtual void insert(libMesh::NumericVector<libMesh::Number> & vector) override; 392 : virtual void insertLower(libMesh::NumericVector<libMesh::Number> & vector) override; 393 : virtual void add(libMesh::NumericVector<libMesh::Number> & vector) override; 394 : 395 : virtual void setActiveTags(const std::set<TagID> & vtags) override; 396 : 397 : [[noreturn]] virtual const MooseArray<OutputType> & nodalValueArray() const override; 398 : [[noreturn]] virtual const MooseArray<OutputType> & nodalValueOldArray() const override; 399 : [[noreturn]] virtual const MooseArray<OutputType> & nodalValueOlderArray() const override; 400 : 401 187 : virtual const FieldVariablePhiValue & phi() const override final { return _phi; } 402 0 : virtual const FieldVariablePhiGradient & gradPhi() const override final { return _grad_phi; } 403 : [[noreturn]] virtual const FieldVariablePhiSecond & secondPhi() const override final; 404 : [[noreturn]] const FieldVariablePhiValue & curlPhi() const override final; 405 : [[noreturn]] const FieldVariablePhiDivergence & divPhi() const override final; 406 : 407 0 : virtual const FieldVariablePhiValue & phiFace() const override final { return _phi_face; } 408 0 : virtual const FieldVariablePhiGradient & gradPhiFace() const override final 409 : { 410 0 : return _grad_phi_face; 411 : } 412 : [[noreturn]] virtual const FieldVariablePhiSecond & secondPhiFace() const override final; 413 : 414 0 : virtual const FieldVariablePhiValue & phiFaceNeighbor() const override final 415 : { 416 0 : return _phi_face_neighbor; 417 : } 418 0 : virtual const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final 419 : { 420 0 : return _grad_phi_face_neighbor; 421 : } 422 : [[noreturn]] virtual const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final; 423 : 424 0 : virtual const FieldVariablePhiValue & phiNeighbor() const override final { return _phi_neighbor; } 425 0 : virtual const FieldVariablePhiGradient & gradPhiNeighbor() const override final 426 : { 427 0 : return _grad_phi_neighbor; 428 : } 429 : [[noreturn]] virtual const FieldVariablePhiSecond & secondPhiNeighbor() const override final; 430 : 431 : virtual const FieldVariableValue & vectorTagValue(TagID tag) const override; 432 : virtual const DofValues & vectorTagDofValue(TagID tag) const override; 433 : [[noreturn]] virtual const DofValues & nodalMatrixTagValue(TagID tag) const override; 434 : virtual const FieldVariableValue & matrixTagValue(TagID tag) const override; 435 : 436 : virtual const FieldVariableValue & sln() const override; 437 : virtual const FieldVariableValue & slnOld() const override; 438 : virtual const FieldVariableValue & slnOlder() const override; 439 : virtual const FieldVariableGradient & gradSln() const override; 440 : virtual const FieldVariableGradient & gradSlnOld() const override; 441 : virtual const FieldVariableValue & slnNeighbor() const override; 442 : virtual const FieldVariableValue & slnOldNeighbor() const override; 443 : virtual const FieldVariableGradient & gradSlnNeighbor() const override; 444 : virtual const FieldVariableGradient & gradSlnOldNeighbor() const override; 445 : 446 : [[noreturn]] virtual const ADTemplateVariableSecond<OutputType> & adSecondSln() const override; 447 : [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adUDot() const override; 448 : [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adUDotDot() const override; 449 : [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> & adGradSlnDot() const override; 450 : [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adSlnNeighbor() const override; 451 : [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> & 452 : adGradSlnNeighbor() const override; 453 : [[noreturn]] virtual const ADTemplateVariableSecond<OutputType> & 454 : adSecondSlnNeighbor() const override; 455 : [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adUDotNeighbor() const override; 456 : [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & 457 : adUDotDotNeighbor() const override; 458 : [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> & 459 : adGradSlnNeighborDot() const override; 460 : [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adSln() const override; 461 : [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> & adGradSln() const override; 462 : [[noreturn]] virtual const ADTemplateVariableCurl<OutputType> & adCurlSln() const override; 463 : [[noreturn]] virtual const ADTemplateVariableCurl<OutputType> & 464 : adCurlSlnNeighbor() const override; 465 : 466 : virtual const DofValues & dofValues() const override; 467 : virtual const DofValues & dofValuesOld() const override; 468 : 469 : virtual const DofValues & dofValuesOlder() const override; 470 : virtual const DofValues & dofValuesPreviousNL() const override; 471 : virtual const DofValues & dofValuesNeighbor() const override; 472 : virtual const DofValues & dofValuesOldNeighbor() const override; 473 : virtual const DofValues & dofValuesOlderNeighbor() const override; 474 : virtual const DofValues & dofValuesPreviousNLNeighbor() const override; 475 : [[noreturn]] virtual const DofValues & dofValuesDot() const override; 476 : [[noreturn]] virtual const DofValues & dofValuesDotNeighbor() const override; 477 : [[noreturn]] virtual const DofValues & dofValuesDotOld() const override; 478 : [[noreturn]] virtual const DofValues & dofValuesDotOldNeighbor() const override; 479 : [[noreturn]] virtual const DofValues & dofValuesDotDot() const override; 480 : [[noreturn]] virtual const DofValues & dofValuesDotDotNeighbor() const override; 481 : [[noreturn]] virtual const DofValues & dofValuesDotDotOld() const override; 482 : [[noreturn]] virtual const DofValues & dofValuesDotDotOldNeighbor() const override; 483 : [[noreturn]] virtual const MooseArray<libMesh::Number> & dofValuesDuDotDu() const override; 484 : [[noreturn]] virtual const MooseArray<libMesh::Number> & 485 : dofValuesDuDotDuNeighbor() const override; 486 : [[noreturn]] virtual const MooseArray<libMesh::Number> & dofValuesDuDotDotDu() const override; 487 : [[noreturn]] virtual const MooseArray<libMesh::Number> & 488 : dofValuesDuDotDotDuNeighbor() const override; 489 : 490 : [[noreturn]] virtual const ADDofValues & adDofValues() const override; 491 : [[noreturn]] virtual const ADDofValues & adDofValuesNeighbor() const override; 492 : [[noreturn]] virtual const ADDofValues & adDofValuesDot() const override; 493 : [[noreturn]] virtual const dof_id_type & nodalDofIndex() const override final; 494 : [[noreturn]] virtual const dof_id_type & nodalDofIndexNeighbor() const override final; 495 : 496 0 : virtual std::size_t phiSize() const override final { return _phi.size(); } 497 0 : virtual std::size_t phiFaceSize() const override final { return _phi_face.size(); } 498 0 : virtual std::size_t phiNeighborSize() const override final { return _phi_neighbor.size(); } 499 0 : virtual std::size_t phiFaceNeighborSize() const override final 500 : { 501 0 : return _phi_face_neighbor.size(); 502 : } 503 : [[noreturn]] virtual std::size_t phiLowerSize() const override final; 504 : }; 505 : 506 : template <typename OutputType> 507 : typename MooseLinearVariableFV<OutputType>::ValueType 508 6858 : MooseLinearVariableFV<OutputType>::evaluate(const ElemArg & elem_arg, const StateArg & state) const 509 : { 510 6858 : const auto & elem_info = this->_mesh.elemInfo(elem_arg.elem->id()); 511 6858 : return getElemValue(elem_info, state); 512 : } 513 : 514 : template <typename OutputType> 515 : typename MooseLinearVariableFV<OutputType>::ValueType 516 0 : MooseLinearVariableFV<OutputType>::evaluate(const ElemPointArg & elem_point, 517 : const StateArg & state) const 518 : { 519 0 : const auto & elem_info = this->_mesh.elemInfo(elem_point.elem->id()); 520 0 : return getElemValue(elem_info, state); 521 : } 522 : 523 : template <typename OutputType> 524 : typename MooseLinearVariableFV<OutputType>::ValueType 525 99552 : MooseLinearVariableFV<OutputType>::evaluate(const ElemQpArg & elem_qp, const StateArg & state) const 526 : { 527 99552 : const auto & elem_info = this->_mesh.elemInfo(elem_qp.elem->id()); 528 99552 : return getElemValue(elem_info, state); 529 : } 530 : 531 : template <typename OutputType> 532 : typename MooseLinearVariableFV<OutputType>::ValueType 533 0 : MooseLinearVariableFV<OutputType>::evaluate(const ElemSideQpArg & elem_side_qp, 534 : const StateArg & state) const 535 : { 536 0 : return (*this)(ElemPointArg{elem_side_qp.elem, elem_side_qp.point, false}, state); 537 : } 538 : 539 : template <typename OutputType> 540 : typename MooseLinearVariableFV<OutputType>::GradientType 541 1400 : MooseLinearVariableFV<OutputType>::evaluateGradient(const ElemQpArg & qp_arg, 542 : const StateArg & state) const 543 : { 544 1400 : const auto & elem_info = this->_mesh.elemInfo(qp_arg.elem->id()); 545 1400 : return gradSln(elem_info, state); 546 : } 547 : 548 : template <typename OutputType> 549 : typename MooseLinearVariableFV<OutputType>::GradientType 550 0 : MooseLinearVariableFV<OutputType>::evaluateGradient(const ElemArg & elem_arg, 551 : const StateArg & state) const 552 : { 553 0 : const auto & elem_info = this->_mesh.elemInfo(elem_arg.elem->id()); 554 0 : return gradSln(elem_info, state); 555 : } 556 : 557 : template <typename OutputType> 558 : typename MooseLinearVariableFV<OutputType>::GradientType 559 140 : MooseLinearVariableFV<OutputType>::evaluateGradient(const FaceArg & face, 560 : const StateArg & state) const 561 : { 562 : mooseAssert(face.fi, "We must have a non-null face information"); 563 140 : return gradSln(*face.fi, state); 564 : } 565 : 566 : template <typename OutputType> 567 : void 568 0 : MooseLinearVariableFV<OutputType>::timeIntegratorError() const 569 : { 570 0 : mooseError("MooseLinearVariableFV does not support time integration at the moment! The variable " 571 : "which is causing the issue: ", 572 0 : this->name()); 573 : } 574 : 575 : template <typename OutputType> 576 : void 577 0 : MooseLinearVariableFV<OutputType>::gradientStateError(const StateArg & state) const 578 : { 579 0 : mooseError("MooseLinearVariableFV does not currently support ElemInfo/FaceInfo gradient " 580 : "evaluation for non-current states. Requested state index ", 581 0 : state.state, 582 : " for variable '", 583 0 : this->name(), 584 : "'. Old-state requests typically use state index 1."); 585 : } 586 : 587 : template <typename OutputType> 588 : void 589 0 : MooseLinearVariableFV<OutputType>::lowerDError() const 590 : { 591 0 : mooseError("Lower dimensional element support not implemented for finite volume variables!The " 592 : "variable which is causing the issue: ", 593 0 : this->name()); 594 : } 595 : 596 : template <typename OutputType> 597 : void 598 0 : MooseLinearVariableFV<OutputType>::nodalError() const 599 : { 600 0 : mooseError("FV variables don't support nodal variable treatment! The variable which is causing " 601 : "the issue: ", 602 0 : this->name()); 603 : } 604 : 605 : template <typename OutputType> 606 : void 607 0 : MooseLinearVariableFV<OutputType>::adError() const 608 : { 609 0 : mooseError("Linear FV variable does not support automatic differentiation, the variable which is " 610 : "attempting it is: ", 611 0 : this->name()); 612 : } 613 : 614 : // Declare all the specializations, as the template specialization declarations below must know 615 : template <> 616 : ADReal MooseLinearVariableFV<Real>::evaluateDot(const ElemArg & elem, const StateArg & state) const; 617 : 618 : // Prevent implicit instantiation in other translation units where these classes are used 619 : extern template class MooseLinearVariableFV<Real>;