https://mooseframework.inl.gov
ADRayKernel.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 // Local Includes
13 #include "IntegralRayKernelBase.h"
14 
15 // MOOSE Includes
16 #include "MooseVariableInterface.h"
17 #include "TaggingInterface.h"
18 
19 // Forward declarations
20 template <typename>
22 
24 // Not implementing this until there is a use case and tests for it!
25 // using ADVectorRayKernel = ADRayKernelTempl<RealVectorValue>;
26 
30 template <typename T>
32  public MooseVariableInterface<T>,
33  public TaggingInterface
34 {
35 public:
36  ADRayKernelTempl(const InputParameters & params);
37 
39 
44 
45  void onSegment() override final;
46 
47 protected:
51  virtual ADReal computeQpResidual() = 0;
52 
56  virtual void precalculateResidual(){};
57 
60 
63 
66 
69 
72 
74  unsigned int _i;
75 
77  unsigned int _j;
78 
79 private:
83  void computeJacobian();
87  void computeResidual();
88 
90  std::vector<ADReal> _residuals;
91 };
std::vector< ADReal > _residuals
Temporary for filling the residuals.
Definition: ADRayKernel.h:90
Base class for a RayKernel that integrates along a Ray segment.
Base class for an AD ray kernel that contributes to the residual and/or Jacobian. ...
Definition: ADRayKernel.h:21
typename OutputTools< T >::VariableTestValue ADTemplateVariableTestValue
unsigned int _i
Current index for the test function.
Definition: ADRayKernel.h:74
void onSegment() override final
Called on each segment of a Ray.
Definition: ADRayKernel.C:68
typename OutputTools< T >::VariablePhiValue ADTemplateVariablePhiValue
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
static InputParameters validParams()
Definition: ADRayKernel.C:22
const ADTemplateVariableValue< T > & _u
Holds the solution at current quadrature points.
Definition: ADRayKernel.h:65
DualNumber< Real, DNDerivativeType, true > ADReal
MooseVariableField< T > & variable()
The MooseVariable this RayKernel contributes to.
Definition: ADRayKernel.h:43
ADRayKernelTempl(const InputParameters &params)
Definition: ADRayKernel.C:34
const ADTemplateVariablePhiValue< T > & _phi
The current shape functions.
Definition: ADRayKernel.h:71
unsigned int _j
Current index for the shape function.
Definition: ADRayKernel.h:77
void computeResidual()
Computes and contributes to the residual for a segment.
Definition: ADRayKernel.C:82
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
MooseVariableField< T > & _var
The MooseVariable this kernel contributes to.
Definition: ADRayKernel.h:56
virtual ADReal computeQpResidual()=0
Compute this kernel&#39;s contribution to the residual at _qp and _i.
const ADTemplateVariableTestValue< T > & _test
The current test function.
Definition: ADRayKernel.h:62
virtual void precalculateResidual()
Insertion point for calculation before the residual computation.
Definition: ADRayKernel.h:56
void computeJacobian()
Computes and contributes to the Jacobian for a segment.
Definition: ADRayKernel.C:96
const ADTemplateVariableGradient< T > & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: ADRayKernel.h:68