Line data Source code
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 "KokkosNodalKernelBase.h"
13 :
14 : namespace Moose
15 : {
16 : namespace Kokkos
17 : {
18 :
19 : /**
20 : * The base class for a user to derive their own Kokkos nodal kernels.
21 : *
22 : * The user should define computeQpResidual(), computeQpJacobian(), and computeQpOffDiagJacobian()
23 : * as inlined public methods in their derived class (not virtual override). The signature of
24 : * computeQpResidual() expected to be defined in the derived class is as follows:
25 : *
26 : * @param qp The dummy quadrature point index (= 0)
27 : * @param datum The AssemblyDatum object of the current thread
28 : * @returns The residual contribution
29 : *
30 : * KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const;
31 : *
32 : * The signatures of computeQpJacobian() and computeOffDiagQpJacobian() can be found in the code
33 : * below, and their definition in the derived class is optional. If they are defined in the derived
34 : * class, they will hide the default definitions in the base class.
35 : */
36 : class NodalKernel : public NodalKernelBase
37 : {
38 : public:
39 : static InputParameters validParams();
40 :
41 : /**
42 : * Constructor
43 : */
44 : NodalKernel(const InputParameters & parameters);
45 :
46 : /**
47 : * Dispatch residual calculation
48 : */
49 : virtual void computeResidual() override;
50 : /**
51 : * Dispatch diagonal and off-diagonal Jacobian calculation
52 : */
53 : virtual void computeJacobian() override;
54 :
55 : /**
56 : * Default methods to prevent compile errors even when these methods were not defined in the
57 : * derived class
58 : */
59 : ///@{
60 : /**
61 : * Compute diagonal Jacobian contribution on a node
62 : * @param qp The dummy quadrature point index (= 0)
63 : * @param datum The AssemblyDatum object of the current thread
64 : * @returns The diagonal Jacobian contribution
65 : */
66 0 : KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* qp */,
67 : AssemblyDatum & /* datum */) const
68 : {
69 0 : return 0;
70 : }
71 : /**
72 : * Compute off-diagonal Jacobian contribution on a node
73 : * @param jvar The variable number for column
74 : * @param qp The dummy quadrature point index (= 0)
75 : * @param datum The AssemblyDatum object of the current thread
76 : * @returns The off-diagonal Jacobian contribution
77 : */
78 0 : KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* jvar */,
79 : const unsigned int /* qp */,
80 : AssemblyDatum & /* datum */) const
81 : {
82 0 : return 0;
83 : }
84 : /**
85 : * Get the function pointer of the default computeQpJacobian()
86 : * @returns The function pointer
87 : */
88 262202 : static auto defaultJacobian() { return &NodalKernel::computeQpJacobian; }
89 : /**
90 : * Get the function pointer of the default computeQpOffDiagJacobian()
91 : * @returns The function pointer
92 : */
93 262202 : static auto defaultOffDiagJacobian() { return &NodalKernel::computeQpOffDiagJacobian; }
94 : ///@}
95 :
96 : /**
97 : * The parallel computation entry functions called by Kokkos
98 : */
99 : ///@{
100 : template <typename Derived>
101 : KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived & kernel) const;
102 : template <typename Derived>
103 : KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid, const Derived & kernel) const;
104 : template <typename Derived>
105 : KOKKOS_FUNCTION void
106 : operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & kernel) const;
107 : ///@}
108 :
109 : protected:
110 : /**
111 : * Current solution at nodes
112 : */
113 : const VariableValue _u;
114 :
115 : private:
116 : /**
117 : * Flag whether this kernel is boundary-restricted
118 : */
119 : const bool _boundary_restricted;
120 : };
121 :
122 : template <typename Derived>
123 : KOKKOS_FUNCTION void
124 1265482 : NodalKernel::operator()(ResidualLoop, const ThreadID tid, const Derived & kernel) const
125 : {
126 1265482 : auto node = _boundary_restricted ? kokkosBoundaryNodeID(tid) : kokkosBlockNodeID(tid);
127 1265482 : auto & sys = kokkosSystem(_kokkos_var.sys());
128 :
129 1265482 : if (!sys.isNodalDefined(node, _kokkos_var.var()))
130 0 : return;
131 :
132 1265482 : AssemblyDatum datum(node, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
133 :
134 1265482 : Real local_re = kernel.computeQpResidual(0, datum);
135 :
136 1265482 : accumulateTaggedNodalResidual(true, local_re, node);
137 : }
138 :
139 : template <typename Derived>
140 : KOKKOS_FUNCTION void
141 194664 : NodalKernel::operator()(JacobianLoop, const ThreadID tid, const Derived & kernel) const
142 : {
143 194664 : auto node = _boundary_restricted ? kokkosBoundaryNodeID(tid) : kokkosBlockNodeID(tid);
144 194664 : auto & sys = kokkosSystem(_kokkos_var.sys());
145 :
146 194664 : if (!sys.isNodalDefined(node, _kokkos_var.var()))
147 0 : return;
148 :
149 194664 : AssemblyDatum datum(node, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
150 :
151 194664 : Real local_ke = kernel.computeQpJacobian(0, datum);
152 :
153 194664 : accumulateTaggedNodalMatrix(true, local_ke, node, _kokkos_var.var());
154 : }
155 :
156 : template <typename Derived>
157 : KOKKOS_FUNCTION void
158 249100 : NodalKernel::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & kernel) const
159 : {
160 249100 : auto node = _boundary_restricted ? kokkosBoundaryNodeID(_thread(tid, 1))
161 249100 : : kokkosBlockNodeID(_thread(tid, 1));
162 249100 : auto & sys = kokkosSystem(_kokkos_var.sys());
163 249100 : auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
164 :
165 249100 : if (!sys.isNodalDefined(node, _kokkos_var.var()))
166 0 : return;
167 :
168 249100 : AssemblyDatum datum(node, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
169 :
170 249100 : Real local_ke = kernel.computeQpOffDiagJacobian(jvar, 0, datum);
171 :
172 249100 : accumulateTaggedNodalMatrix(true, local_ke, node, jvar);
173 : }
174 :
175 : } // namespace Kokkos
176 : } // namespace Moose
|