https://mooseframework.inl.gov
KokkosVariableValue.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "KokkosDatum.h"
13 
14 #include "MooseVariableBase.h"
15 
16 namespace Moose
17 {
18 namespace Kokkos
19 {
20 
24 class VariablePhiValue
26 {
27 public:
35  KOKKOS_FUNCTION Real operator()(ResidualDatum & datum, unsigned int i, unsigned int qp) const
36  {
37  auto & elem = datum.elem();
38  auto side = datum.side();
39  auto fe = datum.jfe();
40 
41  return side == libMesh::invalid_uint
42  ? datum.assembly().getPhi(elem.subdomain, elem.type, fe)(i, qp)
43  : datum.assembly().getPhiFace(elem.subdomain, elem.type, fe)(side)(i, qp);
44  }
45 };
46 
48 {
49 public:
57  KOKKOS_FUNCTION Real3 operator()(ResidualDatum & datum, unsigned int i, unsigned int qp) const
58  {
59  auto & elem = datum.elem();
60  auto side = datum.side();
61  auto fe = datum.jfe();
62 
63  return datum.J(qp) *
64  (side == libMesh::invalid_uint
65  ? datum.assembly().getGradPhi(elem.subdomain, elem.type, fe)(i, qp)
66  : datum.assembly().getGradPhiFace(elem.subdomain, elem.type, fe)(side)(i, qp));
67  }
68 };
69 
71 {
72 public:
80  KOKKOS_FUNCTION Real operator()(ResidualDatum & datum, unsigned int i, unsigned int qp) const
81  {
82  auto & elem = datum.elem();
83  auto side = datum.side();
84  auto fe = datum.ife();
85 
86  return side == libMesh::invalid_uint
87  ? datum.assembly().getPhi(elem.subdomain, elem.type, fe)(i, qp)
88  : datum.assembly().getPhiFace(elem.subdomain, elem.type, fe)(side)(i, qp);
89  }
90 };
91 
93 {
94 public:
102  KOKKOS_FUNCTION Real3 operator()(ResidualDatum & datum, unsigned int i, unsigned int qp) const
103  {
104  auto & elem = datum.elem();
105  auto side = datum.side();
106  auto fe = datum.ife();
107 
108  return datum.J(qp) *
109  (side == libMesh::invalid_uint
110  ? datum.assembly().getGradPhi(elem.subdomain, elem.type, fe)(i, qp)
111  : datum.assembly().getGradPhiFace(elem.subdomain, elem.type, fe)(side)(i, qp));
112  }
113 };
115 
119 class VariableValue
121 {
122 public:
128  VariableValue(Variable var, bool nodal = false) : _var(var), _nodal(nodal)
129  {
130  if (nodal && !var.nodal())
131  mooseError("Cannot get nodal values of a non-nodal variable.");
132  }
140  const TagName & tag = Moose::SOLUTION_TAG,
141  bool nodal = false)
142  : _var(var, tag), _nodal(nodal)
143  {
144  if (nodal && !var.isNodal())
145  mooseError("Cannot get nodal values of a non-nodal variable.");
146  }
147 
152  KOKKOS_FUNCTION operator bool() const { return _var.coupled(); }
153 
161  KOKKOS_FUNCTION Real operator()(Datum & datum, unsigned int qp, unsigned int comp = 0) const
162  {
163  if (_var.coupled())
164  {
165  if (_nodal)
166  {
167  KOKKOS_ASSERT(datum.isNodal());
168 
169  auto node = datum.node();
170  auto dof = datum.system(_var.sys(comp)).getNodeLocalDofIndex(node, _var.var(comp));
171 
172  return datum.system(_var.sys(comp)).getVectorDofValue(dof, _var.tag());
173  }
174  else
175  {
176  KOKKOS_ASSERT(!datum.isNodal());
177 
178  auto & elem = datum.elem();
179  auto side = datum.side();
180  auto qp_offset = datum.qpOffset();
181 
182  return side == libMesh::invalid_uint
183  ? datum.system(_var.sys(comp))
184  .getVectorQpValue(elem, qp_offset + qp, _var.var(comp), _var.tag())
185  : datum.system(_var.sys(comp))
186  .getVectorQpValueFace(elem, side, qp, _var.var(comp), _var.tag());
187  }
188  }
189  else
190  return _var.value(comp);
191  }
192 
193 private:
201  const bool _nodal;
202 };
203 
205 {
206 public:
217  VariableGradient(const MooseVariableBase & var, const TagName & tag = Moose::SOLUTION_TAG)
218  : _var(var, tag)
219  {
220  }
221 
226  KOKKOS_FUNCTION operator bool() const { return _var.coupled(); }
227 
235  KOKKOS_FUNCTION Real3 operator()(Datum & datum, unsigned int qp, unsigned int comp = 0) const
236  {
237  if (_var.coupled())
238  {
239  auto & elem = datum.elem();
240  auto side = datum.side();
241  auto qp_offset = datum.qpOffset();
242 
243  return side == libMesh::invalid_uint
244  ? datum.system(_var.sys(comp))
245  .getVectorQpGrad(elem, qp_offset + qp, _var.var(comp), _var.tag())
246  : datum.system(_var.sys(comp))
247  .getVectorQpGradFace(
248  elem, side, datum.J(qp), qp, _var.var(comp), _var.tag());
249  }
250  else
251  return Real3(0);
252  }
253 
254 private:
259 };
261 
262 } // namespace Kokkos
263 } // namespace Moose
KOKKOS_FUNCTION const auto & getPhi(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the shape functions of a FE type for an element type and subdomain.
KOKKOS_FUNCTION TagID tag() const
Get the vector tag ID.
virtual bool isNodal() const
Is this variable nodal.
Variable _var
Coupled Kokkos variable.
KOKKOS_FUNCTION unsigned int sys(unsigned int comp=0) const
Get the system number of a component.
The Kokkos object that holds thread-private data in the parallel operations of any Kokkos object...
Definition: KokkosDatum.h:25
KOKKOS_FUNCTION Real operator()(Datum &datum, unsigned int qp, unsigned int comp=0) const
Get the current variable value.
const unsigned int invalid_uint
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
OutputTools< Real >::VariablePhiValue VariablePhiValue
Definition: MooseTypes.h:320
KOKKOS_FUNCTION const auto & getGradPhiFace(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the gradient of face shape functions of a FE type for an element type and subdomain.
KOKKOS_FUNCTION const auto & getGradPhi(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the gradient of shape functions of a FE type for an element type and subdomain.
const bool _nodal
Flag whether nodal values are requested.
KOKKOS_FUNCTION const auto & getPhiFace(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the face shape functions of a FE type for an element type and subdomain.
KOKKOS_FUNCTION bool coupled() const
Get whether the variable is coupled.
KOKKOS_FUNCTION const Real33 & J(const unsigned int qp)
Get the inverse of Jacobian matrix | dxi/dx deta/dx dzeta/dx | | dxi/dy deta/dy dzeta/dy | | dxi/dz d...
Definition: KokkosDatum.h:125
KOKKOS_FUNCTION unsigned int side() const
Get the side index.
Definition: KokkosDatum.h:90
KOKKOS_FUNCTION Real3 operator()(Datum &datum, unsigned int qp, unsigned int comp=0) const
Get the current variable gradient.
KOKKOS_FUNCTION Real operator()(ResidualDatum &datum, unsigned int i, unsigned int qp) const
Get the current test function.
VariableValue(Variable var, bool nodal=false)
Constructor.
KOKKOS_FUNCTION unsigned int ife() const
Get the variable FE type ID.
Definition: KokkosDatum.h:340
KOKKOS_FUNCTION const Assembly & assembly() const
Get the Kokkos assembly.
Definition: KokkosDatum.h:68
VariableGradient(const MooseVariableBase &var, const TagName &tag=Moose::SOLUTION_TAG)
Constructor.
VariableValue(const MooseVariableBase &var, const TagName &tag=Moose::SOLUTION_TAG, bool nodal=false)
Constructor.
KOKKOS_FUNCTION Real3 operator()(ResidualDatum &datum, unsigned int i, unsigned int qp) const
Get the gradient of the current test function.
KOKKOS_FUNCTION bool nodal() const
Get whether the variable is nodal.
KOKKOS_FUNCTION Real operator()(ResidualDatum &datum, unsigned int i, unsigned int qp) const
Get the current shape function.
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:315
KOKKOS_FUNCTION ContiguousNodeID node() const
Get the contiguous node ID.
Definition: KokkosDatum.h:95
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION Real value(unsigned int comp=0) const
Get the default value of a component.
KOKKOS_FUNCTION bool isNodal() const
Get whether the current datum is on a node.
Definition: KokkosDatum.h:115
The Kokkos variable object that carries the coupled variable and tag information. ...
KOKKOS_FUNCTION const ElementInfo & elem() const
Get the element information object.
Definition: KokkosDatum.h:80
const TagName SOLUTION_TAG
Definition: MooseTypes.C:25
Variable _var
Coupled Kokkos variable.
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
KOKKOS_FUNCTION unsigned int var(unsigned int comp=0) const
Get the variable number of a component.
VariableGradient(Variable var)
Constructor.
KOKKOS_FUNCTION const System & system(unsigned int sys) const
Get the Kokkos system.
Definition: KokkosDatum.h:74
KOKKOS_FUNCTION dof_id_type qpOffset() const
Get the starting offset into the global quadrature point index.
Definition: KokkosDatum.h:105
KOKKOS_FUNCTION unsigned int jfe() const
Get the coupled variable FE type ID.
Definition: KokkosDatum.h:345
The Kokkos object that holds thread-private data in the parallel operations of Kokkos residual object...
Definition: KokkosDatum.h:245
KOKKOS_FUNCTION Real3 operator()(ResidualDatum &datum, unsigned int i, unsigned int qp) const
Get the gradient of the current shape function.
Base variable class.