https://mooseframework.inl.gov
KokkosDirichletBCBase.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 "KokkosNodalBC.h"
13 
14 namespace Moose
15 {
16 namespace Kokkos
17 {
18 
22 template <typename Derived>
23 class DirichletBCBase : public NodalBC<Derived>
24 {
26 
27 public:
29 
34 
39  virtual bool preset() const override { return _preset; }
40 
45  virtual void presetSolution(TagID tag) override;
46 
50  KOKKOS_FUNCTION void operator()(const ThreadID tid) const;
51 
57  KOKKOS_FUNCTION Real computeQpResidual(const ContiguousNodeID node) const;
58 
59 private:
63  const bool _preset;
68 };
69 
70 template <typename Derived>
73 {
75  params.addParam<bool>(
76  "preset", true, "Whether or not to preset the BC (apply the value before the solve begins).");
77  return params;
78 }
79 
80 template <typename Derived>
82  : NodalBC<Derived>(parameters), _preset(this->template getParam<bool>("preset"))
83 {
84 }
85 
86 template <typename Derived>
87 void
89 {
90  _solution_tag = tag;
91 
92  ::Kokkos::parallel_for(::Kokkos::RangePolicy<ExecSpace, ::Kokkos::IndexType<ThreadID>>(
93  0, this->numKokkosBoundaryNodes()),
94  *static_cast<Derived *>(this));
95 }
96 
97 template <typename Derived>
98 KOKKOS_FUNCTION void
100 {
101  auto bc = static_cast<const Derived *>(this);
102  auto node = kokkosBoundaryNodeID(tid);
103  auto & sys = kokkosSystem(_kokkos_var.sys());
104  auto dof = sys.getNodeLocalDofIndex(node, _kokkos_var.var());
105 
107  return;
108 
109  sys.getVectorDofValue(dof, _solution_tag) = bc->computeValue(node);
110 }
111 
112 template <typename Derived>
113 KOKKOS_FUNCTION Real
115 {
116  auto bc = static_cast<const Derived *>(this);
117 
118  return _u(node) - bc->computeValue(node);
119 }
120 
121 } // namespace Kokkos
122 } // namespace Moose
123 
124 #define usingKokkosDirichletBCBaseMembers(T) \
125  usingKokkosNodalBCMembers(T); \
126  \
127 public: \
128  using Moose::Kokkos::DirichletBCBase<T>::operator()
virtual bool preset() const override
Get whether the value is to be preset.
const bool _preset
Flag whether the value is to be preset.
unsigned int TagID
Definition: MooseTypes.h:210
dof_id_type ThreadID
Definition: KokkosThread.h:18
virtual void presetSolution(TagID tag) override
Dispatch solution vector preset.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
virtual Real computeQpResidual() override
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
TagID _solution_tag
Tag associated with the solution vector to be preset.
DirichletBCBase(const InputParameters &parameters)
Constructor.
static InputParameters validParams()
dof_id_type ContiguousNodeID
Definition: KokkosMesh.h:21
static const dof_id_type invalid_id
static InputParameters validParams()
KOKKOS_FUNCTION void operator()(const ThreadID tid) const
The preset function called by Kokkos.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
KOKKOS_FUNCTION Real computeQpResidual(const ContiguousNodeID node) const
Compute residual contribution on a node.
The base class for a user to derive their own Kokkos nodal boundary conditions.
Definition: KokkosNodalBC.h:49
The base Kokkos boundary condition of a Dirichlet type.