Line data Source code
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 :
19 : /**
20 : * The base class for a user to derive their own Kokkos integrated boundary conditions.
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 i The element-local DOF index
27 : * @param qp The local quadrature point index
28 : * @param datum The AssemblyDatum object of the current thread
29 : * @returns The residual contribution
30 : *
31 : * KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i,
32 : * const unsigned int qp,
33 : * AssemblyDatum & datum) const;
34 : *
35 : * The signatures of computeQpJacobian() and computeOffDiagQpJacobian() can be found in the code
36 : * below, and their definition in the derived class is optional. If they are defined in the derived
37 : * class, they will hide the default definitions in the base class.
38 : */
39 : class IntegratedBC : public IntegratedBCBase
40 : {
41 : public:
42 : static InputParameters validParams();
43 :
44 : /**
45 : * Constructor
46 : */
47 : IntegratedBC(const InputParameters & parameters);
48 :
49 : /**
50 : * Dispatch residual calculation
51 : */
52 : virtual void computeResidual() override;
53 : /**
54 : * Dispatch diagonal and off-diagonal Jacobian calculation
55 : */
56 : virtual void computeJacobian() override;
57 :
58 : /**
59 : * Default methods to prevent compile errors even when these methods were not defined in the
60 : * derived class
61 : */
62 : ///@{
63 : /**
64 : * Compute diagonal Jacobian contribution on a quadrature point
65 : * @param i The test function DOF index
66 : * @param j The trial function DOF index
67 : * @param qp The local quadrature point index
68 : * @param datum The AssemblyDatum object of the current thread
69 : * @returns The diagonal Jacobian contribution
70 : */
71 0 : KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
72 : const unsigned int /* j */,
73 : const unsigned int /* qp */,
74 : AssemblyDatum & /* datum */) const
75 : {
76 0 : return 0;
77 : }
78 : /**
79 : * Compute off-diagonal Jacobian contribution on a quadrature point
80 : * @param i The test function DOF index
81 : * @param j The trial function DOF index
82 : * @param jvar The variable number for column
83 : * @param qp The local quadrature point index
84 : * @param datum The AssemblyDatum object of the current thread
85 : * @returns The off-diagonal Jacobian contribution
86 : */
87 0 : 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 0 : return 0;
94 : }
95 : /**
96 : * Get the function pointer of the default computeQpJacobian()
97 : * @returns The function pointer
98 : */
99 225928 : static auto defaultJacobian() { return &IntegratedBC::computeQpJacobian; }
100 : /**
101 : * Get the function pointer of the default computeQpOffDiagJacobian()
102 : * @returns The function pointer
103 : */
104 225928 : static auto defaultOffDiagJacobian() { return &IntegratedBC::computeQpOffDiagJacobian; }
105 : ///@}
106 :
107 : /**
108 : * Shims for hook methods that can be leveraged to implement static polymorphism
109 : */
110 : ///{@
111 : template <typename Derived>
112 845484 : KOKKOS_FUNCTION Real computeQpResidualShim(const Derived & bc,
113 : const unsigned int i,
114 : const unsigned int qp,
115 : AssemblyDatum & datum) const
116 : {
117 845484 : return bc.computeQpResidual(i, qp, datum);
118 : }
119 : template <typename Derived>
120 42880 : 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 42880 : return bc.computeQpJacobian(i, j, qp, datum);
127 : }
128 : template <typename Derived>
129 80 : 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 80 : return bc.computeQpOffDiagJacobian(i, j, jvar, qp, datum);
137 : }
138 : ///@}
139 :
140 : /**
141 : * The parallel computation entry functions called by Kokkos
142 : */
143 : ///@{
144 : 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;
151 : ///@}
152 :
153 : /**
154 : * The parallel computation bodies that can be customized in the derived class by defining
155 : * them in the derived class with the same signature.
156 : * Make sure to define them as inlined public methods if to be defined in the derived class.
157 : */
158 : ///@{
159 : /**
160 : * Compute residual
161 : * @param bc The boundary condition object of the final derived type
162 : * @param datum The AssemblyDatum object of the current thread
163 : */
164 : template <typename Derived>
165 : KOKKOS_FUNCTION void computeResidualInternal(const Derived & bc, AssemblyDatum & datum) const;
166 : /**
167 : * Compute diagonal Jacobian
168 : * @param bc The boundary condition object of the final derived type
169 : * @param datum The AssemblyDatum object of the current thread
170 : */
171 : template <typename Derived>
172 : KOKKOS_FUNCTION void computeJacobianInternal(const Derived & bc, AssemblyDatum & datum) const;
173 : /**
174 : * Compute off-diagonal Jacobian
175 : * @param bc The boundary condition object of the final derived type
176 : * @param datum The AssemblyDatum object of the current thread
177 : */
178 : template <typename Derived>
179 : KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived & bc,
180 : AssemblyDatum & datum) const;
181 : ///@}
182 :
183 : protected:
184 : /**
185 : * Current test function
186 : */
187 : const VariableTestValue _test;
188 : /**
189 : * Gradient of the current test function
190 : */
191 : const VariableTestGradient _grad_test;
192 : /**
193 : * Current shape function
194 : */
195 : const VariablePhiValue _phi;
196 : /**
197 : * Gradient of the current shape function
198 : */
199 : const VariablePhiGradient _grad_phi;
200 : /**
201 : * Current solution at quadrature points
202 : */
203 : const VariableValue _u;
204 : /**
205 : * Gradient of the current solution at quadrature points
206 : */
207 : const VariableGradient _grad_u;
208 : };
209 :
210 : template <typename Derived>
211 : KOKKOS_FUNCTION void
212 98420 : IntegratedBC::operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const
213 : {
214 98420 : auto [elem, side] = kokkosBoundaryElementSideID(tid);
215 :
216 196840 : AssemblyDatum datum(
217 98420 : elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
218 :
219 98420 : bc.computeResidualInternal(bc, datum);
220 98420 : }
221 :
222 : template <typename Derived>
223 : KOKKOS_FUNCTION void
224 1340 : IntegratedBC::operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const
225 : {
226 1340 : auto [elem, side] = kokkosBoundaryElementSideID(tid);
227 :
228 2680 : AssemblyDatum datum(
229 1340 : elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
230 :
231 1340 : bc.computeJacobianInternal(bc, datum);
232 1340 : }
233 :
234 : template <typename Derived>
235 : KOKKOS_FUNCTION void
236 20 : IntegratedBC::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const
237 : {
238 20 : auto [elem, side] = kokkosBoundaryElementSideID(_thread(tid, 1));
239 :
240 20 : auto & sys = kokkosSystem(_kokkos_var.sys());
241 20 : auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
242 :
243 20 : if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
244 0 : return;
245 :
246 20 : AssemblyDatum datum(elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
247 :
248 20 : bc.computeOffDiagJacobianInternal(bc, datum);
249 : }
250 :
251 : template <typename Derived>
252 : KOKKOS_FUNCTION void
253 98420 : IntegratedBC::computeResidualInternal(const Derived & bc, AssemblyDatum & datum) const
254 : {
255 98420 : ResidualObject::computeResidualInternal(
256 : datum,
257 196840 : [&](Real * local_re, const unsigned int ib, const unsigned int ie)
258 : {
259 310044 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
260 : {
261 211624 : datum.reinit();
262 :
263 1057108 : for (unsigned int i = ib; i < ie; ++i)
264 845484 : local_re[i] += datum.JxW(qp) * bc.computeQpResidualShim(bc, i, qp, datum);
265 : }
266 : });
267 98420 : }
268 :
269 : template <typename Derived>
270 : KOKKOS_FUNCTION void
271 1340 : IntegratedBC::computeJacobianInternal(const Derived & bc, AssemblyDatum & datum) const
272 : {
273 1340 : ResidualObject::computeJacobianInternal(
274 : datum,
275 2680 : [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
276 : {
277 4020 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
278 : {
279 2680 : datum.reinit();
280 :
281 45560 : for (unsigned int ij = ijb; ij < ije; ++ij)
282 : {
283 42880 : unsigned int i = ij % datum.n_jdofs();
284 42880 : unsigned int j = ij / datum.n_jdofs();
285 :
286 42880 : local_ke[ij] += datum.JxW(qp) * bc.computeQpJacobianShim(bc, i, j, qp, datum);
287 : }
288 : }
289 : });
290 1340 : }
291 :
292 : template <typename Derived>
293 : KOKKOS_FUNCTION void
294 20 : IntegratedBC::computeOffDiagJacobianInternal(const Derived & bc, AssemblyDatum & datum) const
295 : {
296 20 : ResidualObject::computeJacobianInternal(
297 : datum,
298 40 : [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
299 : {
300 40 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
301 : {
302 20 : datum.reinit();
303 :
304 100 : for (unsigned int ij = ijb; ij < ije; ++ij)
305 : {
306 80 : unsigned int i = ij % datum.n_jdofs();
307 80 : unsigned int j = ij / datum.n_jdofs();
308 :
309 80 : local_ke[ij] +=
310 80 : datum.JxW(qp) * bc.computeQpOffDiagJacobianShim(bc, i, j, datum.jvar(), qp, datum);
311 : }
312 : }
313 : });
314 20 : }
315 :
316 : } // namespace Kokkos
317 : } // namespace Moose
|