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 
52 template <typename Derived>
54 {
55 public:
57 
62 
66  virtual void computeResidual() override;
70  virtual void computeJacobian() override;
71 
76 
85  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
86  const unsigned int /* j */,
87  const unsigned int /* qp */,
88  ResidualDatum & /* datum */) const
89  {
90  return 0;
91  }
101  KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* i */,
102  const unsigned int /* j */,
103  const unsigned int /* jvar */,
104  const unsigned int /* qp */,
105  ResidualDatum & /* datum */) const
106  {
107  return 0;
108  }
110 
114  KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid) const;
116  KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid) const;
117  KOKKOS_FUNCTION void operator()(OffDiagJacobianLoop, const ThreadID tid) const;
119 
125 
131  KOKKOS_FUNCTION void computeResidualInternal(const Derived * bc, ResidualDatum & datum) const;
137  KOKKOS_FUNCTION void computeJacobianInternal(const Derived * bc, ResidualDatum & datum) const;
143  KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived * bc,
144  ResidualDatum & datum) const;
146 
147 protected:
172 
173 protected:
178  virtual bool defaultJacobian() const
179  {
180  return &Derived::computeQpJacobian == &IntegratedBC::computeQpJacobian;
181  }
186  virtual bool defaultOffDiagJacobian() const
187  {
188  return &Derived::computeQpOffDiagJacobian == &IntegratedBC::computeQpOffDiagJacobian;
189  }
190 };
191 
192 template <typename Derived>
195 {
197  return params;
198 }
199 
200 template <typename Derived>
203  _test(),
204  _grad_test(),
205  _phi(),
206  _grad_phi(),
207  _u(_var),
208  _grad_u(_var)
209 {
211 }
212 
213 template <typename Derived>
214 void
216 {
217  ::Kokkos::RangePolicy<ResidualLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
218  0, numKokkosBoundarySides());
219  ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
220  ::Kokkos::fence();
221 }
222 
223 template <typename Derived>
224 void
226 {
227  if (!defaultJacobian())
228  {
229  ::Kokkos::RangePolicy<JacobianLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
230  0, numKokkosBoundarySides());
231  ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
232  ::Kokkos::fence();
233  }
234 
235  if (!defaultOffDiagJacobian())
236  {
237  auto & sys = kokkosSystem(_kokkos_var.sys());
238 
239  _thread.resize({sys.getCoupling(_kokkos_var.var()).size(), numKokkosBoundarySides()});
240 
241  ::Kokkos::RangePolicy<OffDiagJacobianLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
242  0, _thread.size());
243  ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
244  ::Kokkos::fence();
245  }
246 }
247 
248 template <typename Derived>
249 KOKKOS_FUNCTION void
251 {
252  auto bc = static_cast<const Derived *>(this);
253  auto [elem, side] = kokkosBoundaryElementSideID(tid);
254 
255  ResidualDatum datum(
256  elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
257 
258  bc->computeResidualInternal(bc, datum);
259 }
260 
261 template <typename Derived>
262 KOKKOS_FUNCTION void
264 {
265  auto bc = static_cast<const Derived *>(this);
266  auto [elem, side] = kokkosBoundaryElementSideID(tid);
267 
268  ResidualDatum datum(
269  elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
270 
271  bc->computeJacobianInternal(bc, datum);
272 }
273 
274 template <typename Derived>
275 KOKKOS_FUNCTION void
277 {
278  auto bc = static_cast<const Derived *>(this);
279  auto [elem, side] = kokkosBoundaryElementSideID(_thread(tid, 1));
280 
281  auto & sys = kokkosSystem(_kokkos_var.sys());
282  auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
283 
284  if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
285  return;
286 
287  ResidualDatum datum(elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
288 
289  bc->computeOffDiagJacobianInternal(bc, datum);
290 }
291 
292 template <typename Derived>
293 KOKKOS_FUNCTION void
295 {
296  ResidualObject::computeResidualInternal(
297  datum,
298  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
299  {
300  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
301  {
302  datum.reinit();
303 
304  for (unsigned int i = ib; i < ie; ++i)
305  local_re[i] += datum.JxW(qp) * bc->computeQpResidual(i, qp, datum);
306  }
307  });
308 }
309 
310 template <typename Derived>
311 KOKKOS_FUNCTION void
313 {
314  ResidualObject::computeJacobianInternal(
315  datum,
316  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
317  {
318  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
319  {
320  datum.reinit();
321 
322  for (unsigned int ij = ijb; ij < ije; ++ij)
323  {
324  unsigned int i = ij % datum.n_jdofs();
325  unsigned int j = ij / datum.n_jdofs();
326 
327  local_ke[ij] += datum.JxW(qp) * bc->computeQpJacobian(i, j, qp, datum);
328  }
329  }
330  });
331 }
332 
333 template <typename Derived>
334 KOKKOS_FUNCTION void
336  ResidualDatum & datum) const
337 {
338  ResidualObject::computeJacobianInternal(
339  datum,
340  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
341  {
342  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
343  {
344  datum.reinit();
345 
346  for (unsigned int ij = ijb; ij < ije; ++ij)
347  {
348  unsigned int i = ij % datum.n_jdofs();
349  unsigned int j = ij / datum.n_jdofs();
350 
351  local_ke[ij] +=
352  datum.JxW(qp) * bc->computeQpOffDiagJacobian(i, j, datum.jvar(), qp, datum);
353  }
354  }
355  });
356 }
357 
358 } // namespace Kokkos
359 } // namespace Moose
360 
361 #define usingKokkosIntegratedBCMembers(T) \
362  usingKokkosIntegratedBCBaseMembers; \
363  \
364 protected: \
365  using Moose::Kokkos::IntegratedBC<T>::_test; \
366  using Moose::Kokkos::IntegratedBC<T>::_grad_test; \
367  using Moose::Kokkos::IntegratedBC<T>::_phi; \
368  using Moose::Kokkos::IntegratedBC<T>::_grad_phi; \
369  using Moose::Kokkos::IntegratedBC<T>::_u; \
370  using Moose::Kokkos::IntegratedBC<T>::_grad_u; \
371  \
372 public: \
373  using Moose::Kokkos::IntegratedBC<T>::operator()
VarFieldType
Definition: MooseTypes.h:722
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:287
const VariableGradient _grad_u
Gradient of the current solution at quadrature points.
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int, const unsigned int, const unsigned int, const unsigned int, ResidualDatum &) const
Compute off-diagonal Jacobian contribution on a quadrature point.
dof_id_type ThreadID
Definition: KokkosThread.h:18
The base class for a user to derive their own Kokkos integrated boundary conditions.
KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived *bc, ResidualDatum &datum) const
Compute off-diagonal Jacobian.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, const unsigned int, const unsigned int, ResidualDatum &) const
Default methods to prevent compile errors even when these methods were not defined in the derived cla...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The base class for Kokkos integrated boundary conditions.
IntegratedBC(const InputParameters &parameters)
Constructor.
const VariableTestValue _test
Current test function.
The Kokkos wrapper classes for MOOSE-like shape function access.
const VariablePhiValue _phi
Current shape function.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: IntegratedBC.C:109
virtual bool defaultJacobian() const
Get whether computeQpJacobian() was not defined in the derived class.
virtual void computeResidual() override
Dispatch residual calculation.
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:95
static InputParameters validParams()
KOKKOS_FUNCTION void computeResidualInternal(const Derived *bc, ResidualDatum &datum) const
The parallel computation bodies that can be customized in the derived class by defining them in the d...
MooseVariable & _var
Definition: IntegratedBC.h:82
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
const VariableTestGradient _grad_test
Gradient of the current test function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION void computeJacobianInternal(const Derived *bc, ResidualDatum &datum) const
Compute diagonal Jacobian.
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: IntegratedBC.C:85
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:126
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
static InputParameters validParams()
const VariablePhiGradient _grad_phi
Gradient of the current shape function.
virtual bool defaultOffDiagJacobian() const
Get whether computeQpOffDiagJacobian() was not defined in the derived class.
The Kokkos wrapper classes for MOOSE-like variable value access.
const VariableValue _u
Current solution at quadrature points.
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
Definition: KokkosDatum.h:302
The Kokkos object that holds thread-private data in the parallel operations of Kokkos residual object...
Definition: KokkosDatum.h:222
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Definition: KokkosDatum.h:147
KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid) const
The parallel computation entry functions called by Kokkos.