https://mooseframework.inl.gov
KokkosIntegratedBC.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 "KokkosIntegratedBCBase.h"
13 
14 namespace Moose
15 {
16 namespace Kokkos
17 {
18 
40 {
41 public:
43 
48 
52  virtual void computeResidual() override;
56  virtual void computeJacobian() override;
57 
62 
71  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
72  const unsigned int /* j */,
73  const unsigned int /* qp */,
74  AssemblyDatum & /* datum */) const
75  {
76  return 0;
77  }
87  KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* i */,
88  const unsigned int /* j */,
89  const unsigned int /* jvar */,
90  const unsigned int /* qp */,
91  AssemblyDatum & /* datum */) const
92  {
93  return 0;
94  }
106 
110  template <typename Derived>
112  KOKKOS_FUNCTION Real computeQpResidualShim(const Derived & bc,
113  const unsigned int i,
114  const unsigned int qp,
115  AssemblyDatum & datum) const
116  {
117  return bc.computeQpResidual(i, qp, datum);
118  }
119  template <typename Derived>
120  KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived & bc,
121  const unsigned int i,
122  const unsigned int j,
123  const unsigned int qp,
124  AssemblyDatum & datum) const
125  {
126  return bc.computeQpJacobian(i, j, qp, datum);
127  }
128  template <typename Derived>
129  KOKKOS_FUNCTION Real computeQpOffDiagJacobianShim(const Derived & bc,
130  const unsigned int i,
131  const unsigned int j,
132  const unsigned int jvar,
133  const unsigned int qp,
134  AssemblyDatum & datum) const
135  {
136  return bc.computeQpOffDiagJacobian(i, j, jvar, qp, datum);
137  }
139 
143  template <typename Derived>
145  KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const;
146  template <typename Derived>
147  KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const;
148  template <typename Derived>
149  KOKKOS_FUNCTION void
150  operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const;
152 
158 
164  template <typename Derived>
165  KOKKOS_FUNCTION void computeResidualInternal(const Derived & bc, AssemblyDatum & datum) const;
171  template <typename Derived>
172  KOKKOS_FUNCTION void computeJacobianInternal(const Derived & bc, AssemblyDatum & datum) const;
178  template <typename Derived>
179  KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived & bc,
180  AssemblyDatum & datum) const;
182 
183 protected:
208 };
209 
210 template <typename Derived>
211 KOKKOS_FUNCTION void
212 IntegratedBC::operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const
213 {
214  auto [elem, side] = kokkosBoundaryElementSideID(tid);
215 
216  AssemblyDatum datum(
218 
219  bc.computeResidualInternal(bc, datum);
220 }
221 
222 template <typename Derived>
223 KOKKOS_FUNCTION void
224 IntegratedBC::operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const
225 {
226  auto [elem, side] = kokkosBoundaryElementSideID(tid);
227 
228  AssemblyDatum datum(
230 
231  bc.computeJacobianInternal(bc, datum);
232 }
233 
234 template <typename Derived>
235 KOKKOS_FUNCTION void
236 IntegratedBC::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const
237 {
238  auto [elem, side] = kokkosBoundaryElementSideID(_thread(tid, 1));
239 
240  auto & sys = kokkosSystem(_kokkos_var.sys());
241  auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
242 
243  if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
244  return;
245 
246  AssemblyDatum datum(elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
247 
248  bc.computeOffDiagJacobianInternal(bc, datum);
249 }
250 
251 template <typename Derived>
252 KOKKOS_FUNCTION void
253 IntegratedBC::computeResidualInternal(const Derived & bc, AssemblyDatum & datum) const
254 {
256  datum,
257  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
258  {
259  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
260  {
261  datum.reinit();
262 
263  for (unsigned int i = ib; i < ie; ++i)
264  local_re[i] += datum.JxW(qp) * bc.computeQpResidualShim(bc, i, qp, datum);
265  }
266  });
267 }
268 
269 template <typename Derived>
270 KOKKOS_FUNCTION void
271 IntegratedBC::computeJacobianInternal(const Derived & bc, AssemblyDatum & datum) const
272 {
274  datum,
275  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
276  {
277  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
278  {
279  datum.reinit();
280 
281  for (unsigned int ij = ijb; ij < ije; ++ij)
282  {
283  unsigned int i = ij % datum.n_jdofs();
284  unsigned int j = ij / datum.n_jdofs();
285 
286  local_ke[ij] += datum.JxW(qp) * bc.computeQpJacobianShim(bc, i, j, qp, datum);
287  }
288  }
289  });
290 }
291 
292 template <typename Derived>
293 KOKKOS_FUNCTION void
295 {
297  datum,
298  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
299  {
300  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
301  {
302  datum.reinit();
303 
304  for (unsigned int ij = ijb; ij < ije; ++ij)
305  {
306  unsigned int i = ij % datum.n_jdofs();
307  unsigned int j = ij / datum.n_jdofs();
308 
309  local_ke[ij] +=
310  datum.JxW(qp) * bc.computeQpOffDiagJacobianShim(bc, i, j, datum.jvar(), qp, datum);
311  }
312  }
313  });
314 }
315 
316 } // namespace Kokkos
317 } // namespace Moose
KOKKOS_FUNCTION unsigned int sys(unsigned int comp=0) const
Get the system number of a component.
KOKKOS_FUNCTION void computeJacobianInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental Jacobian.
KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived &bc, const unsigned int i, const unsigned int j, const unsigned int qp, AssemblyDatum &datum) const
KOKKOS_FUNCTION void computeResidualInternal(const Derived &bc, AssemblyDatum &datum) const
The parallel computation bodies that can be customized in the derived class by defining them in the d...
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
static auto defaultJacobian()
Get the function pointer of the default computeQpJacobian()
dof_id_type ThreadID
Definition: KokkosThread.h:18
The base class for a user to derive their own Kokkos integrated boundary conditions.
KOKKOS_FUNCTION Real computeQpResidualShim(const Derived &bc, const unsigned int i, const unsigned int qp, AssemblyDatum &datum) const
Shims for hook methods that can be leveraged to implement static polymorphism.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Thread _thread
Kokkos thread object.
The base class for Kokkos integrated boundary conditions.
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
Definition: KokkosDatum.h:329
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:314
const VariableValue _u
Current solution at quadrature points.
The Kokkos wrapper classes for MOOSE-like shape function access.
const VariableGradient _grad_u
Gradient of the current solution at quadrature points.
const VariablePhiGradient _grad_phi
Gradient of the current shape function.
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int, const unsigned int, const unsigned int, const unsigned int, AssemblyDatum &) const
Compute off-diagonal Jacobian contribution on a quadrature point.
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:100
KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived &bc) const
The parallel computation entry functions called by Kokkos.
KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
Get the const reference of a Kokkos system.
Definition: KokkosSystem.h:625
const VariableTestGradient _grad_test
Gradient of the current test function.
static auto defaultOffDiagJacobian()
Get the function pointer of the default computeQpOffDiagJacobian()
const VariableTestValue _test
Current test function.
KOKKOS_FUNCTION const Mesh & kokkosMesh() const
Get the const reference of the Kokkos mesh.
Definition: KokkosMesh.h:360
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived &bc, AssemblyDatum &datum) const
Compute off-diagonal Jacobian.
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &bc, AssemblyDatum &datum) const
Compute diagonal Jacobian.
KOKKOS_FUNCTION Real computeQpOffDiagJacobianShim(const Derived &bc, const unsigned int i, const unsigned int j, const unsigned int jvar, const unsigned int qp, AssemblyDatum &datum) const
KOKKOS_FUNCTION auto kokkosBoundaryElementSideID(ThreadID tid) const
Get the contiguous element ID - side index pair this Kokkos thread is operating on.
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
const VariablePhiValue _phi
Current shape function.
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
Definition: KokkosDatum.h:244
virtual void computeResidual() override
Dispatch residual calculation.
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:139
Variable _kokkos_var
Kokkos variable.
static InputParameters validParams()
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
The Kokkos wrapper classes for MOOSE-like variable value access.
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, const unsigned int, const unsigned int, AssemblyDatum &) const
Default methods to prevent compile errors even when these methods were not defined in the derived cla...
KOKKOS_FUNCTION unsigned int var(unsigned int comp=0) const
Get the variable number of a component.
KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const
Get the list of off-diagonal coupled variable numbers of a variable.
Definition: KokkosSystem.h:147
KOKKOS_FUNCTION void computeResidualInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental residual.
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Definition: KokkosDatum.h:166
IntegratedBC(const InputParameters &parameters)
Constructor.
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.
Definition: KokkosSystem.h:608