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 "ADKernel.h" 13 : 14 : /** 15 : * This ADKernel adds standardized methods for assembling to a primary 16 : * scalar variable associated with the primary variable of the ADKernel 17 : * object. Essentially, the entire row of the residual and Jacobian 18 : * associated with this scalar variable will also be assembled here 19 : * using the loops over volumetric elements. 20 : * This variable is "scalar_variable" in the input file and "kappa" 21 : * within the source code. 22 : */ 23 : class ADKernelScalarBase : public ADKernel 24 : { 25 : public: 26 : static InputParameters validParams(); 27 : 28 : ADKernelScalarBase(const InputParameters & parameters); 29 : 30 : /** 31 : * The scalar variable that this kernel operates on. 32 : */ 33 : const MooseVariableScalar & scalarVariable() const 34 : { 35 : mooseAssert(_kappa_var_ptr, "kappa pointer should have been set in the constructor"); 36 : return *_kappa_var_ptr; 37 : } 38 : 39 : virtual void computeResidual() override; 40 : virtual void computeJacobian() override; 41 : /** 42 : * Computes d-_var-residual / d-jvar as well as d-_kappa-residual / d-jvar 43 : */ 44 : virtual void computeOffDiagJacobian(unsigned int jvar) override; 45 : /** 46 : * Computes jacobian block with respect to a scalar variable 47 : * @param jvar, the number of the (other) scalar variable 48 : */ 49 : void computeOffDiagJacobianScalar(unsigned int /*jvar_num*/) override; 50 : 51 : /** 52 : * Computes residual and jacobian block for field and scalar variables 53 : */ 54 : void computeResidualAndJacobian() override; 55 : 56 : protected: 57 : /** 58 : * Method for computing the scalar part of residual at quadrature points 59 : */ 60 : virtual ADReal computeScalarQpResidual(); 61 : 62 : /** 63 : * compute the \p _scalar_residuals member for filling the Jacobian. We want to calculate these 64 : * residuals up-front when doing loal derivative indexing because we can use those residuals to 65 : * fill \p _local_ke for every associated jvariable. We do not want to re-do these calculations 66 : * for every jvariable and corresponding \p _local_ke. For global indexing we will simply pass 67 : * the computed \p _residuals directly to \p Assembly::addJacobian 68 : */ 69 : virtual void computeScalarResidualsForJacobian(); 70 : 71 : /** 72 : * For coupling scalar variables 73 : * Added solely for GenericKernelScalar override; should not be used 74 : */ 75 0 : virtual Real computeQpOffDiagJacobianScalar(unsigned int /*jvar*/) { return 0; } 76 : 77 : /** 78 : * Method for computing the scalar variable part of Jacobian at quadrature points 79 : * Added solely for GenericKernelScalar override; should not be used 80 : */ 81 0 : virtual Real computeScalarQpJacobian() { return 0; } 82 : 83 : /** 84 : * Method for computing an off-diagonal jacobian component at quadrature points. 85 : * Added solely for GenericKernelScalar override; should not be used 86 : */ 87 0 : virtual Real computeScalarQpOffDiagJacobian(const unsigned int /*jvar_num*/) { return 0; } 88 : 89 : /** 90 : * Method for computing an off-diagonal jacobian component at quadrature points. 91 : * Added solely for GenericKernelScalar override; should not be used 92 : */ 93 0 : virtual Real computeScalarQpOffDiagJacobianScalar(const unsigned int /*svar_num*/) { return 0; } 94 : 95 : /** 96 : * Put necessary evaluations depending on qp but independent of test functions here 97 : */ 98 16275 : virtual void initScalarQpResidual() {} 99 : 100 : /// Whether a scalar variable is declared for this kernel 101 : const bool _use_scalar; 102 : 103 : /// Whether to compute scalar contributions for this instance 104 : const bool _compute_scalar_residuals; 105 : 106 : /// Whether to compute field contributions for this instance 107 : const bool _compute_field_residuals; 108 : 109 : /// (Pointer to) Scalar variable this kernel operates on 110 : const MooseVariableScalar * const _kappa_var_ptr; 111 : 112 : /// The unknown scalar variable ID 113 : const unsigned int _kappa_var; 114 : 115 : /// Order of the scalar variable, used in several places 116 : const unsigned int _k_order; 117 : 118 : /// Reference to the current solution at the current quadrature point 119 : const ADVariableValue & _kappa; 120 : 121 : /// Used internally to iterate over each scalar component 122 : unsigned int _h; 123 : unsigned int _l; 124 : std::vector<ADReal> _scalar_residuals; 125 : }; 126 : 127 : inline ADReal 128 3 : ADKernelScalarBase::computeScalarQpResidual() 129 : { 130 3 : mooseError( 131 : "A scalar_variable has been set and compute_scalar_residuals=true, ", 132 : "but the computeScalarQpResidual method was not overridden. Accidental call of base class?"); 133 : return 0; 134 : }