https://mooseframework.inl.gov
RayKernel.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 "RayKernel.h"
11 
12 // Local includes
13 #include "RayTracingStudy.h"
14 
15 // MOOSE includes
16 #include "Assembly.h"
17 #include "NonlinearSystemBase.h"
18 
19 template <typename T>
22 {
23  auto params = IntegralRayKernelBase::validParams();
25 
26  params.template addRequiredParam<NonlinearVariableName>(
27  "variable", "The name of the variable that this RayKernel operates on");
28 
29  return params;
30 }
31 
32 template <typename T>
34  : IntegralRayKernelBase(params),
35  MooseVariableInterface<T>(this,
36  false,
37  "variable",
41  TaggingInterface(this),
42  _var(*this->mooseVariable()),
43  _u(_is_implicit ? _var.sln() : _var.slnOld()),
44  _grad_u(_is_implicit ? _var.gradSln() : _var.gradSlnOld()),
45  _test(_var.phi()),
46  _grad_test(_var.gradPhi()),
47  _phi(_assembly.phi(_var)),
48  _grad_phi(_assembly.gradPhi(_var))
49 {
50  // We do not allow RZ/RSPHERICAL because in the context of these coord
51  // systems there is no way to represent a line source - we would end up
52  // with a plane/surface source or a volumetric source, respectively.
53  // This is also why we do not multiply by _coord[_qp] in any of the
54  // integrations that follow.
55  for (const auto & subdomain_id : _mesh.meshSubdomains())
56  if (_fe_problem.getCoordSystem(subdomain_id) != Moose::COORD_XYZ)
57  mooseError("Not valid on coordinate systems other than XYZ");
58 
60 }
61 
62 template <typename T>
63 void
65 {
66  mooseAssert(_current_subdomain_id == _assembly.currentSubdomainID(), "Subdomain IDs not in sync");
67  mooseAssert(_fe_problem.currentlyComputingJacobian() || _fe_problem.currentlyComputingResidual(),
68  "Not computing residual or Jacobian");
69 
70  if (_fe_problem.currentlyComputingJacobian())
71  computeJacobian();
72  else if (_fe_problem.currentlyComputingResidual())
73  computeResidual();
74 }
75 
76 template <typename T>
77 void
79 {
80  prepareVectorTag(_assembly, _var.number());
81 
82  precalculateResidual();
83  for (_i = 0; _i < _test.size(); ++_i)
84  for (_qp = 0; _qp < _q_point.size(); ++_qp)
85  _local_re(_i) += _JxW[_qp] * computeQpResidual();
86 
87  accumulateTaggedLocalResidual();
88 }
89 
90 template <typename T>
91 void
93 {
94  if (!isImplicit())
95  return;
96 
97  _subproblem.prepareShapes(_var.number(), _tid);
98 
99  precalculateJacobian();
100 
101  const auto & ce = _fe_problem.couplingEntries(_tid, _nl->number());
102  for (const auto & it : ce)
103  {
104  MooseVariableFEBase & ivariable = *(it.first);
105  MooseVariableFEBase & jvariable = *(it.second);
106 
107  const auto ivar = ivariable.number();
108  if (ivar != _var.number() || !jvariable.activeOnSubdomain(_current_subdomain_id))
109  continue;
110 
111  // This error should be caught by the coupleable interface
112  mooseAssert(jvariable.count() == 1, "ArrayMooseVariable objects are not coupleable");
113 
114  const auto jvar = jvariable.number();
115 
116  prepareMatrixTag(_assembly, _var.number(), jvar);
117  if (_local_ke.m() != _test.size())
118  return;
119 
120  if (ivar == jvar)
121  {
122  precalculateJacobian();
123  for (_i = 0; _i < _test.size(); _i++)
124  for (_j = 0; _j < _phi.size(); _j++)
125  for (_qp = 0; _qp < _JxW.size(); _qp++)
126  _local_ke(_i, _j) += _JxW[_qp] * computeQpJacobian();
127  }
128  else
129  {
130  precalculateOffDiagJacobian(jvar);
131  for (_i = 0; _i < _test.size(); _i++)
132  for (_j = 0; _j < _phi.size(); _j++)
133  for (_qp = 0; _qp < _JxW.size(); _qp++)
134  _local_ke(_i, _j) += _JxW[_qp] * computeQpOffDiagJacobian(jvar);
135  }
136 
137  accumulateTaggedLocalMatrix();
138  }
139 }
140 
141 template class RayKernelTempl<Real>;
142 // Not implementing this until there is a use case and tests for it!
143 // template class RayKernelTempl<RealVectorValue>;
VarFieldType
Base class for a RayKernel that integrates along a Ray segment.
VAR_SOLVER
unsigned int number() const
static InputParameters validParams()
Definition: RayKernel.C:21
unsigned int count() const
FEProblemBase & _fe_problem
The FEProblemBase.
MooseMesh & _mesh
The MooseMesh.
Base class for a ray kernel that contributes to the residual and/or Jacobian.
Definition: RayKernel.h:21
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
VarKindType
static InputParameters validParams()
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
void addMooseVariableDependency(MooseVariableFieldBase *var)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
RayKernelTempl(const InputParameters &params)
Definition: RayKernel.C:33
void mooseError(Args &&... args) const
Moose::CoordinateSystemType getCoordSystem(SubdomainID sid) const
bool activeOnSubdomain(SubdomainID subdomain) const
VAR_FIELD_STANDARD
static InputParameters validParams()
VAR_FIELD_VECTOR
void computeResidual()
Computes and contributes to the residual for a segment.
Definition: RayKernel.C:78
void computeJacobian()
Computes and contributes to the Jacobian for a segment.
Definition: RayKernel.C:92
const std::set< SubdomainID > & meshSubdomains() const