https://mooseframework.inl.gov
RayKernel.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 VectorRayKernel = RayKernelTempl<RealVectorValue>;
26 
30 template <typename T>
32  public MooseVariableInterface<T>,
33  public TaggingInterface
34 {
35 public:
36  RayKernelTempl(const InputParameters & params);
37 
39 
44 
45  void onSegment() override final;
46 
47 protected:
51  virtual Real computeQpResidual() = 0;
55  virtual Real computeQpJacobian() { return 0; }
60  virtual Real computeQpOffDiagJacobian(const unsigned int /* jvar_num */) { return 0; }
61 
65  virtual void precalculateResidual(){};
69  virtual void precalculateJacobian(){};
73  virtual void precalculateOffDiagJacobian(unsigned int /* jvar_num */){};
74 
77 
80 
83 
86 
89 
92 
95 
97  unsigned int _i;
98 
100  unsigned int _j;
101 
102 private:
106  void computeJacobian();
110  void computeResidual();
111 };
112 
113 extern template class RayKernelTempl<Real>;
virtual Real computeQpJacobian()
Compute this RayKernel&#39;s contribution to the residual at _qp, _i, and _j.
Definition: RayKernel.h:55
Base class for a RayKernel that integrates along a Ray segment.
virtual Real computeQpResidual()=0
Compute this RayKernel&#39;s contribution to the residual at _qp and _i.
static InputParameters validParams()
Definition: RayKernel.C:21
virtual void precalculateJacobian()
Insertion point for calculation before the Jacobian contribution.
Definition: RayKernel.h:69
const OutputTools< T >::VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: RayKernel.h:82
MooseVariableFE< T > & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: RayKernel.h:73
Base class for a ray kernel that contributes to the residual and/or Jacobian.
Definition: RayKernel.h:21
const OutputTools< T >::VariablePhiValue & _phi
Current shape functions.
Definition: RayKernel.h:91
unsigned int _i
Current index for the test function.
Definition: RayKernel.h:97
void onSegment() override final
Called on each segment of a Ray.
Definition: RayKernel.C:64
MooseVariableFE< T > & variable()
The MooseVariable this RayKernel contributes to.
Definition: RayKernel.h:43
const OutputTools< T >::VariableTestValue & _test
Test function.
Definition: RayKernel.h:85
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
RayKernelTempl(const InputParameters &params)
Definition: RayKernel.C:33
const OutputTools< T >::VariableTestGradient & _grad_test
Gradient of the test function.
Definition: RayKernel.h:88
unsigned int _j
Current index for the shape function.
Definition: RayKernel.h:100
const OutputTools< T >::VariableValue & _u
Holds the solution at current quadrature points.
Definition: RayKernel.h:79
virtual void precalculateResidual()
Insertion point for calculation before the residual contribution.
Definition: RayKernel.h:65
virtual void precalculateOffDiagJacobian(unsigned int)
Insertion point for calculation before an off-diagonal Jacobian contribution.
Definition: RayKernel.h:73
void computeResidual()
Computes and contributes to the residual for a segment.
Definition: RayKernel.C:78
const OutputTools< T >::VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: RayKernel.h:94
void computeJacobian()
Computes and contributes to the Jacobian for a segment.
Definition: RayKernel.C:92
virtual Real computeQpOffDiagJacobian(const unsigned int)
Compute this RayKernel&#39;s contribution to the off diagonal jacobian for the variable numbered jvar_num...
Definition: RayKernel.h:60