https://mooseframework.inl.gov
ADScalarKernel.C
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 #include "ADScalarKernel.h"
11 
12 #include "Assembly.h"
13 
16 {
18  return params;
19 }
20 
22  : ScalarKernelBase(parameters), _u(_var.adSln()), _jacobian_already_computed(false)
23 {
24 }
25 
26 void
28 {
30 }
31 
32 void
34 {
36 
37  for (_i = 0; _i < _var.order(); _i++)
39 
41 }
42 
43 void
45 {
47 }
48 
49 void
51 {
53  {
56  }
57 }
58 
59 void
61 {
62 }
63 
64 void
66 {
67  if (_residuals.size() != _var.order())
68  _residuals.resize(_var.order(), 0.0);
69  for (_i = 0; _i < _var.order(); _i++)
71 
73 }
ADScalarKernel(const InputParameters &parameters)
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
std::vector< ADReal > _residuals
Residuals for each order.
unsigned int number() const
Get variable number coming from libMesh.
bool _jacobian_already_computed
Flag indicating that the Jacobian has already been computed.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
MooseVariableScalar & _var
Scalar variable.
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided residual derivatives into the Jacobian for the provided dof indices.
void computeADJacobian()
Computes the Jacobian using automatic differentiation.
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
static InputParameters validParams()
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
libMesh::Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
virtual void reinit() override
Reinitialization method called before each call to computeResidual()
virtual ADReal computeQpResidual()=0
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
Base class shared by AD and non-AD scalar kernels.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
static InputParameters validParams()