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 224716 : static auto defaultJacobian() { return &IntegratedBC::computeQpJacobian; }
100 : /**
101 : * Get the function pointer of the default computeQpOffDiagJacobian()
102 : * @returns The function pointer
103 : */
104 224716 : static auto defaultOffDiagJacobian() { return &IntegratedBC::computeQpOffDiagJacobian; }
105 : ///@}
106 :
107 : /**
108 : * The parallel computation entry functions called by Kokkos
109 : */
110 : ///@{
111 : template <typename Derived>
112 : KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const;
113 : template <typename Derived>
114 : KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const;
115 : template <typename Derived>
116 : KOKKOS_FUNCTION void
117 : operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const;
118 : ///@}
119 :
120 : /**
121 : * The parallel computation bodies that can be customized in the derived class by defining
122 : * them in the derived class with the same signature.
123 : * Make sure to define them as inlined public methods if to be defined in the derived class.
124 : */
125 : ///@{
126 : /**
127 : * Compute residual
128 : * @param bc The boundary condition object of the final derived type
129 : * @param datum The AssemblyDatum object of the current thread
130 : */
131 : template <typename Derived>
132 : KOKKOS_FUNCTION void computeResidualInternal(const Derived & bc, AssemblyDatum & datum) const;
133 : /**
134 : * Compute diagonal Jacobian
135 : * @param bc The boundary condition object of the final derived type
136 : * @param datum The AssemblyDatum object of the current thread
137 : */
138 : template <typename Derived>
139 : KOKKOS_FUNCTION void computeJacobianInternal(const Derived & bc, AssemblyDatum & datum) const;
140 : /**
141 : * Compute off-diagonal Jacobian
142 : * @param bc The boundary condition object of the final derived type
143 : * @param datum The AssemblyDatum object of the current thread
144 : */
145 : template <typename Derived>
146 : KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived & bc,
147 : AssemblyDatum & datum) const;
148 : ///@}
149 :
150 : protected:
151 : /**
152 : * Current test function
153 : */
154 : const VariableTestValue _test;
155 : /**
156 : * Gradient of the current test function
157 : */
158 : const VariableTestGradient _grad_test;
159 : /**
160 : * Current shape function
161 : */
162 : const VariablePhiValue _phi;
163 : /**
164 : * Gradient of the current shape function
165 : */
166 : const VariablePhiGradient _grad_phi;
167 : /**
168 : * Current solution at quadrature points
169 : */
170 : const VariableValue _u;
171 : /**
172 : * Gradient of the current solution at quadrature points
173 : */
174 : const VariableGradient _grad_u;
175 : };
176 :
177 : template <typename Derived>
178 : KOKKOS_FUNCTION void
179 43250 : IntegratedBC::operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const
180 : {
181 43250 : auto [elem, side] = kokkosBoundaryElementSideID(tid);
182 :
183 86500 : AssemblyDatum datum(
184 43250 : elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
185 :
186 43250 : bc.computeResidualInternal(bc, datum);
187 43250 : }
188 :
189 : template <typename Derived>
190 : KOKKOS_FUNCTION void
191 1340 : IntegratedBC::operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const
192 : {
193 1340 : auto [elem, side] = kokkosBoundaryElementSideID(tid);
194 :
195 2680 : AssemblyDatum datum(
196 1340 : elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
197 :
198 1340 : bc.computeJacobianInternal(bc, datum);
199 1340 : }
200 :
201 : template <typename Derived>
202 : KOKKOS_FUNCTION void
203 20 : IntegratedBC::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const
204 : {
205 20 : auto [elem, side] = kokkosBoundaryElementSideID(_thread(tid, 1));
206 :
207 20 : auto & sys = kokkosSystem(_kokkos_var.sys());
208 20 : auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
209 :
210 20 : if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
211 0 : return;
212 :
213 20 : AssemblyDatum datum(elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
214 :
215 20 : bc.computeOffDiagJacobianInternal(bc, datum);
216 : }
217 :
218 : template <typename Derived>
219 : KOKKOS_FUNCTION void
220 43250 : IntegratedBC::computeResidualInternal(const Derived & bc, AssemblyDatum & datum) const
221 : {
222 43250 : ResidualObject::computeResidualInternal(
223 : datum,
224 86500 : [&](Real * local_re, const unsigned int ib, const unsigned int ie)
225 : {
226 144534 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
227 : {
228 101284 : datum.reinit();
229 :
230 505408 : for (unsigned int i = ib; i < ie; ++i)
231 404124 : local_re[i] += datum.JxW(qp) * bc.computeQpResidual(i, qp, datum);
232 : }
233 : });
234 43250 : }
235 :
236 : template <typename Derived>
237 : KOKKOS_FUNCTION void
238 1340 : IntegratedBC::computeJacobianInternal(const Derived & bc, AssemblyDatum & datum) const
239 : {
240 1340 : ResidualObject::computeJacobianInternal(
241 : datum,
242 2680 : [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
243 : {
244 4020 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
245 : {
246 2680 : datum.reinit();
247 :
248 45560 : for (unsigned int ij = ijb; ij < ije; ++ij)
249 : {
250 42880 : unsigned int i = ij % datum.n_jdofs();
251 42880 : unsigned int j = ij / datum.n_jdofs();
252 :
253 42880 : local_ke[ij] += datum.JxW(qp) * bc.computeQpJacobian(i, j, qp, datum);
254 : }
255 : }
256 : });
257 1340 : }
258 :
259 : template <typename Derived>
260 : KOKKOS_FUNCTION void
261 20 : IntegratedBC::computeOffDiagJacobianInternal(const Derived & bc, AssemblyDatum & datum) const
262 : {
263 20 : ResidualObject::computeJacobianInternal(
264 : datum,
265 40 : [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
266 : {
267 40 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
268 : {
269 20 : datum.reinit();
270 :
271 100 : for (unsigned int ij = ijb; ij < ije; ++ij)
272 : {
273 80 : unsigned int i = ij % datum.n_jdofs();
274 80 : unsigned int j = ij / datum.n_jdofs();
275 :
276 80 : local_ke[ij] +=
277 80 : datum.JxW(qp) * bc.computeQpOffDiagJacobian(i, j, datum.jvar(), qp, datum);
278 : }
279 : }
280 : });
281 20 : }
282 :
283 : } // namespace Kokkos
284 : } // namespace Moose
|