https://mooseframework.inl.gov
KokkosVacuumBC.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 #include "KokkosIntegratedBC.h"
13 
20 class KokkosVacuumBC final : public Moose::Kokkos::IntegratedBC<KokkosVacuumBC>
21 {
22 public:
24 
26 
27  KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i,
28  const unsigned int qp,
29  ResidualDatum & datum) const;
30  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i,
31  const unsigned int j,
32  const unsigned int qp,
33  ResidualDatum & datum) const;
34 
35 private:
37  const Real _alpha;
38 };
39 
40 KOKKOS_FUNCTION inline Real
41 KokkosVacuumBC::computeQpResidual(const unsigned int i,
42  const unsigned int qp,
43  ResidualDatum & datum) const
44 {
45  return _test(datum, i, qp) * _alpha * _u(datum, qp) / 2.;
46 }
47 
48 KOKKOS_FUNCTION inline Real
49 KokkosVacuumBC::computeQpJacobian(const unsigned int i,
50  const unsigned int j,
51  const unsigned int qp,
52  ResidualDatum & datum) const
53 {
54  return _test(datum, i, qp) * _alpha * _phi(datum, j, qp) / 2.;
55 }
KokkosVacuumBC(const InputParameters &parameters)
The base class for a user to derive their own Kokkos integrated boundary conditions.
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i, const unsigned int j, const unsigned int qp, ResidualDatum &datum) const
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i, const unsigned int qp, ResidualDatum &datum) const
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...
const VariableTestValue _test
Current test function.
const VariablePhiValue _phi
Current shape function.
const Real _alpha
Ratio of u to du/dn.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
const VariableValue _u
Current solution at quadrature points.
The Kokkos object that holds thread-private data in the parallel operations of Kokkos residual object...
Definition: KokkosDatum.h:222
Implements a simple Vacuum BC for neutron diffusion on the boundary.