1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED: http://github.com/ceed
7
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <assert.h>
12 #include <cuda.h>
13 #include <cuda_runtime.h>
14 #include <stdbool.h>
15 #include <string.h>
16
17 #include "../cuda/ceed-cuda-common.h"
18 #include "../cuda/ceed-cuda-compile.h"
19 #include "ceed-cuda-ref.h"
20
21 //------------------------------------------------------------------------------
22 // Destroy operator
23 //------------------------------------------------------------------------------
CeedOperatorDestroy_Cuda(CeedOperator op)24 static int CeedOperatorDestroy_Cuda(CeedOperator op) {
25 CeedOperator_Cuda *impl;
26
27 CeedCallBackend(CeedOperatorGetData(op, &impl));
28
29 // Apply data
30 CeedCallBackend(CeedFree(&impl->num_points));
31 CeedCallBackend(CeedFree(&impl->skip_rstr_in));
32 CeedCallBackend(CeedFree(&impl->skip_rstr_out));
33 CeedCallBackend(CeedFree(&impl->apply_add_basis_out));
34 CeedCallBackend(CeedFree(&impl->input_field_order));
35 CeedCallBackend(CeedFree(&impl->output_field_order));
36 CeedCallBackend(CeedFree(&impl->input_states));
37
38 for (CeedInt i = 0; i < impl->num_inputs; i++) {
39 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_in[i]));
40 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
41 }
42 CeedCallBackend(CeedFree(&impl->e_vecs_in));
43 CeedCallBackend(CeedFree(&impl->q_vecs_in));
44
45 for (CeedInt i = 0; i < impl->num_outputs; i++) {
46 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_out[i]));
47 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
48 }
49 CeedCallBackend(CeedFree(&impl->e_vecs_out));
50 CeedCallBackend(CeedFree(&impl->q_vecs_out));
51 CeedCallBackend(CeedVectorDestroy(&impl->point_coords_elem));
52
53 // QFunction assembly data
54 for (CeedInt i = 0; i < impl->num_active_in; i++) {
55 CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
56 }
57 CeedCallBackend(CeedFree(&impl->qf_active_in));
58
59 // Diag data
60 if (impl->diag) {
61 Ceed ceed;
62
63 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
64 if (impl->diag->module) {
65 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
66 }
67 if (impl->diag->module_point_block) {
68 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block));
69 }
70 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in));
71 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out));
72 CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
73 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in));
74 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out));
75 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in));
76 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out));
77 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in));
78 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out));
79 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in));
80 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out));
81 CeedCallBackend(CeedDestroy(&ceed));
82 CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
83 CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
84 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr));
85 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
86 }
87 CeedCallBackend(CeedFree(&impl->diag));
88
89 if (impl->asmb) {
90 Ceed ceed;
91
92 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
93 CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
94 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
95 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
96 CeedCallBackend(CeedDestroy(&ceed));
97 }
98 CeedCallBackend(CeedFree(&impl->asmb));
99
100 CeedCallBackend(CeedFree(&impl));
101 return CEED_ERROR_SUCCESS;
102 }
103
104 //------------------------------------------------------------------------------
105 // Setup infields or outfields
106 //------------------------------------------------------------------------------
CeedOperatorSetupFields_Cuda(CeedQFunction qf,CeedOperator op,bool is_input,bool is_at_points,bool * skip_rstr,bool * apply_add_basis,CeedVector * e_vecs,CeedVector * q_vecs,CeedInt num_fields,CeedInt Q,CeedInt num_elem)107 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, bool is_at_points, bool *skip_rstr, bool *apply_add_basis,
108 CeedVector *e_vecs, CeedVector *q_vecs, CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
109 Ceed ceed;
110 CeedQFunctionField *qf_fields;
111 CeedOperatorField *op_fields;
112
113 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
114 if (is_input) {
115 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
116 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
117 } else {
118 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
119 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
120 }
121
122 // Loop over fields
123 for (CeedInt i = 0; i < num_fields; i++) {
124 bool is_active = false, is_strided = false, skip_e_vec = false;
125 CeedSize q_size;
126 CeedInt size;
127 CeedEvalMode eval_mode;
128 CeedVector l_vec;
129 CeedElemRestriction elem_rstr;
130
131 // Check whether this field can skip the element restriction:
132 // Input CEED_VECTOR_ACTIVE
133 // Output CEED_VECTOR_ACTIVE without CEED_EVAL_NONE
134 // Input CEED_VECTOR_NONE with CEED_EVAL_WEIGHT
135 // Input passive vector with CEED_EVAL_NONE and strided restriction with CEED_STRIDES_BACKEND
136 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &l_vec));
137 is_active = l_vec == CEED_VECTOR_ACTIVE;
138 CeedCallBackend(CeedVectorDestroy(&l_vec));
139 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr));
140 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
141 skip_e_vec = (is_input && is_active) || (is_active && eval_mode != CEED_EVAL_NONE) || (eval_mode == CEED_EVAL_WEIGHT);
142 if (!skip_e_vec && is_input && !is_active && eval_mode == CEED_EVAL_NONE) {
143 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
144 if (is_strided) CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_e_vec));
145 }
146 if (skip_e_vec) {
147 e_vecs[i] = NULL;
148 } else {
149 CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i]));
150 }
151 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
152
153 switch (eval_mode) {
154 case CEED_EVAL_NONE:
155 case CEED_EVAL_INTERP:
156 case CEED_EVAL_GRAD:
157 case CEED_EVAL_DIV:
158 case CEED_EVAL_CURL:
159 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
160 q_size = (CeedSize)num_elem * (CeedSize)Q * (CeedSize)size;
161 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
162 break;
163 case CEED_EVAL_WEIGHT: {
164 CeedBasis basis;
165
166 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
167 q_size = (CeedSize)num_elem * (CeedSize)Q;
168 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
169 if (is_at_points) {
170 CeedInt num_points[num_elem];
171
172 for (CeedInt i = 0; i < num_elem; i++) num_points[i] = Q;
173 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, CEED_VECTOR_NONE,
174 q_vecs[i]));
175 } else {
176 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
177 }
178 CeedCallBackend(CeedBasisDestroy(&basis));
179 break;
180 }
181 }
182 }
183 // Drop duplicate restrictions
184 if (is_input) {
185 for (CeedInt i = 0; i < num_fields; i++) {
186 CeedVector vec_i;
187 CeedElemRestriction rstr_i;
188
189 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i));
190 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i));
191 for (CeedInt j = i + 1; j < num_fields; j++) {
192 CeedVector vec_j;
193 CeedElemRestriction rstr_j;
194
195 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j));
196 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j));
197 if (vec_i == vec_j && rstr_i == rstr_j) {
198 if (e_vecs[i]) CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j]));
199 skip_rstr[j] = true;
200 }
201 CeedCallBackend(CeedVectorDestroy(&vec_j));
202 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
203 }
204 CeedCallBackend(CeedVectorDestroy(&vec_i));
205 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
206 }
207 } else {
208 for (CeedInt i = num_fields - 1; i >= 0; i--) {
209 CeedVector vec_i;
210 CeedElemRestriction rstr_i;
211
212 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i));
213 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i));
214 for (CeedInt j = i - 1; j >= 0; j--) {
215 CeedVector vec_j;
216 CeedElemRestriction rstr_j;
217
218 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j));
219 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j));
220 if (vec_i == vec_j && rstr_i == rstr_j) {
221 if (e_vecs[i]) CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j]));
222 skip_rstr[j] = true;
223 apply_add_basis[i] = true;
224 }
225 CeedCallBackend(CeedVectorDestroy(&vec_j));
226 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
227 }
228 CeedCallBackend(CeedVectorDestroy(&vec_i));
229 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
230 }
231 }
232 CeedCallBackend(CeedDestroy(&ceed));
233 return CEED_ERROR_SUCCESS;
234 }
235
236 //------------------------------------------------------------------------------
237 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
238 //------------------------------------------------------------------------------
CeedOperatorSetup_Cuda(CeedOperator op)239 static int CeedOperatorSetup_Cuda(CeedOperator op) {
240 bool is_setup_done;
241 CeedInt Q, num_elem, num_input_fields, num_output_fields;
242 CeedQFunctionField *qf_input_fields, *qf_output_fields;
243 CeedQFunction qf;
244 CeedOperatorField *op_input_fields, *op_output_fields;
245 CeedOperator_Cuda *impl;
246
247 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
248 if (is_setup_done) return CEED_ERROR_SUCCESS;
249
250 CeedCallBackend(CeedOperatorGetData(op, &impl));
251 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
252 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
253 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
254 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
255 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
256
257 // Allocate
258 CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in));
259 CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out));
260 CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in));
261 CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out));
262 CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out));
263 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order));
264 CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order));
265 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states));
266 CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in));
267 CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out));
268 impl->num_inputs = num_input_fields;
269 impl->num_outputs = num_output_fields;
270
271 // Set up infield and outfield e-vecs and q-vecs
272 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, false, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields, Q,
273 num_elem));
274 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, false, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out,
275 impl->q_vecs_out, num_output_fields, Q, num_elem));
276
277 // Reorder fields to allow reuse of buffers
278 impl->max_active_e_vec_len = 0;
279 {
280 bool is_ordered[CEED_FIELD_MAX];
281 CeedInt curr_index = 0;
282
283 for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
284 for (CeedInt i = 0; i < num_input_fields; i++) {
285 CeedSize e_vec_len_i;
286 CeedVector vec_i;
287 CeedElemRestriction rstr_i;
288
289 if (is_ordered[i]) continue;
290 is_ordered[i] = true;
291 impl->input_field_order[curr_index] = i;
292 curr_index++;
293 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
294 if (vec_i == CEED_VECTOR_NONE) {
295 // CEED_EVAL_WEIGHT
296 CeedCallBackend(CeedVectorDestroy(&vec_i));
297 continue;
298 };
299 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
300 CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
301 impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
302 for (CeedInt j = i + 1; j < num_input_fields; j++) {
303 CeedVector vec_j;
304 CeedElemRestriction rstr_j;
305
306 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
307 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
308 if (rstr_i == rstr_j && vec_i == vec_j) {
309 is_ordered[j] = true;
310 impl->input_field_order[curr_index] = j;
311 curr_index++;
312 }
313 CeedCallBackend(CeedVectorDestroy(&vec_j));
314 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
315 }
316 CeedCallBackend(CeedVectorDestroy(&vec_i));
317 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
318 }
319 }
320 {
321 bool is_ordered[CEED_FIELD_MAX];
322 CeedInt curr_index = 0;
323
324 for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false;
325 for (CeedInt i = 0; i < num_output_fields; i++) {
326 CeedSize e_vec_len_i;
327 CeedVector vec_i;
328 CeedElemRestriction rstr_i;
329
330 if (is_ordered[i]) continue;
331 is_ordered[i] = true;
332 impl->output_field_order[curr_index] = i;
333 curr_index++;
334 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i));
335 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i));
336 CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
337 impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
338 for (CeedInt j = i + 1; j < num_output_fields; j++) {
339 CeedVector vec_j;
340 CeedElemRestriction rstr_j;
341
342 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j));
343 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j));
344 if (rstr_i == rstr_j && vec_i == vec_j) {
345 is_ordered[j] = true;
346 impl->output_field_order[curr_index] = j;
347 curr_index++;
348 }
349 CeedCallBackend(CeedVectorDestroy(&vec_j));
350 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
351 }
352 CeedCallBackend(CeedVectorDestroy(&vec_i));
353 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
354 }
355 }
356 CeedCallBackend(CeedClearWorkVectors(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len));
357 {
358 // Create two work vectors for diagonal assembly
359 CeedVector temp_1, temp_2;
360
361 CeedCallBackend(CeedGetWorkVector(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len, &temp_1));
362 CeedCallBackend(CeedGetWorkVector(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len, &temp_2));
363 CeedCallBackend(CeedRestoreWorkVector(CeedOperatorReturnCeed(op), &temp_1));
364 CeedCallBackend(CeedRestoreWorkVector(CeedOperatorReturnCeed(op), &temp_2));
365 }
366 CeedCallBackend(CeedOperatorSetSetupDone(op));
367 CeedCallBackend(CeedQFunctionDestroy(&qf));
368 return CEED_ERROR_SUCCESS;
369 }
370
371 //------------------------------------------------------------------------------
372 // Restrict Operator Inputs
373 //------------------------------------------------------------------------------
CeedOperatorInputRestrict_Cuda(CeedOperatorField op_input_field,CeedQFunctionField qf_input_field,CeedInt input_field,CeedVector in_vec,CeedVector active_e_vec,const bool skip_active,CeedOperator_Cuda * impl,CeedRequest * request)374 static inline int CeedOperatorInputRestrict_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
375 CeedVector in_vec, CeedVector active_e_vec, const bool skip_active, CeedOperator_Cuda *impl,
376 CeedRequest *request) {
377 bool is_active = false;
378 CeedVector l_vec, e_vec = impl->e_vecs_in[input_field];
379
380 // Get input vector
381 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
382 is_active = l_vec == CEED_VECTOR_ACTIVE;
383 if (is_active && skip_active) return CEED_ERROR_SUCCESS;
384 if (is_active) {
385 l_vec = in_vec;
386 if (!e_vec) e_vec = active_e_vec;
387 }
388
389 // Restriction action
390 if (e_vec) {
391 // Restrict, if necessary
392 if (!impl->skip_rstr_in[input_field]) {
393 uint64_t state;
394
395 CeedCallBackend(CeedVectorGetState(l_vec, &state));
396 if (is_active || state != impl->input_states[input_field]) {
397 CeedElemRestriction elem_rstr;
398
399 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_field, &elem_rstr));
400 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, l_vec, e_vec, request));
401 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
402 }
403 impl->input_states[input_field] = state;
404 }
405 }
406 if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
407 return CEED_ERROR_SUCCESS;
408 }
409
410 //------------------------------------------------------------------------------
411 // Input Basis Action
412 //------------------------------------------------------------------------------
CeedOperatorInputBasis_Cuda(CeedOperatorField op_input_field,CeedQFunctionField qf_input_field,CeedInt input_field,CeedVector in_vec,CeedVector active_e_vec,CeedInt num_elem,const bool skip_active,CeedOperator_Cuda * impl)413 static inline int CeedOperatorInputBasis_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
414 CeedVector in_vec, CeedVector active_e_vec, CeedInt num_elem, const bool skip_active,
415 CeedOperator_Cuda *impl) {
416 bool is_active = false;
417 CeedEvalMode eval_mode;
418 CeedVector l_vec, e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field];
419
420 // Skip active input
421 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
422 is_active = l_vec == CEED_VECTOR_ACTIVE;
423 if (is_active && skip_active) return CEED_ERROR_SUCCESS;
424 if (is_active) {
425 l_vec = in_vec;
426 if (!e_vec) e_vec = active_e_vec;
427 }
428
429 // Basis action
430 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
431 switch (eval_mode) {
432 case CEED_EVAL_NONE: {
433 const CeedScalar *e_vec_array;
434
435 if (e_vec) {
436 CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, &e_vec_array));
437 } else {
438 CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, &e_vec_array));
439 }
440 CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
441 break;
442 }
443 case CEED_EVAL_INTERP:
444 case CEED_EVAL_GRAD:
445 case CEED_EVAL_DIV:
446 case CEED_EVAL_CURL: {
447 CeedBasis basis;
448
449 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis));
450 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, e_vec, q_vec));
451 CeedCallBackend(CeedBasisDestroy(&basis));
452 break;
453 }
454 case CEED_EVAL_WEIGHT:
455 break; // No action
456 }
457 if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
458 return CEED_ERROR_SUCCESS;
459 }
460
461 //------------------------------------------------------------------------------
462 // Restore Input Vectors
463 //------------------------------------------------------------------------------
CeedOperatorInputRestore_Cuda(CeedOperatorField op_input_field,CeedQFunctionField qf_input_field,CeedInt input_field,CeedVector in_vec,CeedVector active_e_vec,const bool skip_active,CeedOperator_Cuda * impl)464 static inline int CeedOperatorInputRestore_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
465 CeedVector in_vec, CeedVector active_e_vec, const bool skip_active, CeedOperator_Cuda *impl) {
466 bool is_active = false;
467 CeedEvalMode eval_mode;
468 CeedVector l_vec, e_vec = impl->e_vecs_in[input_field];
469
470 // Skip active input
471 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
472 is_active = l_vec == CEED_VECTOR_ACTIVE;
473 if (is_active && skip_active) return CEED_ERROR_SUCCESS;
474 if (is_active) {
475 l_vec = in_vec;
476 if (!e_vec) e_vec = active_e_vec;
477 }
478
479 // Restore e-vec
480 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
481 if (eval_mode == CEED_EVAL_NONE) {
482 const CeedScalar *e_vec_array;
483
484 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_in[input_field], CEED_MEM_DEVICE, (CeedScalar **)&e_vec_array));
485 if (e_vec) {
486 CeedCallBackend(CeedVectorRestoreArrayRead(e_vec, &e_vec_array));
487 } else {
488 CeedCallBackend(CeedVectorRestoreArrayRead(l_vec, &e_vec_array));
489 }
490 }
491 if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
492 return CEED_ERROR_SUCCESS;
493 }
494
495 //------------------------------------------------------------------------------
496 // Apply and add to output
497 //------------------------------------------------------------------------------
CeedOperatorApplyAdd_Cuda(CeedOperator op,CeedVector in_vec,CeedVector out_vec,CeedRequest * request)498 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
499 CeedInt Q, num_elem, num_input_fields, num_output_fields;
500 Ceed ceed;
501 CeedVector active_e_vec;
502 CeedQFunctionField *qf_input_fields, *qf_output_fields;
503 CeedQFunction qf;
504 CeedOperatorField *op_input_fields, *op_output_fields;
505 CeedOperator_Cuda *impl;
506
507 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
508 CeedCallBackend(CeedOperatorGetData(op, &impl));
509 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
510 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
511 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
512 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
513 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
514
515 // Setup
516 CeedCallBackend(CeedOperatorSetup_Cuda(op));
517
518 // Work vector
519 CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec));
520
521 // Process inputs
522 for (CeedInt i = 0; i < num_input_fields; i++) {
523 CeedInt field = impl->input_field_order[i];
524
525 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, false, impl,
526 request));
527 CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, num_elem, false, impl));
528 }
529
530 // Output pointers, as necessary
531 for (CeedInt i = 0; i < num_output_fields; i++) {
532 CeedEvalMode eval_mode;
533
534 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
535 if (eval_mode == CEED_EVAL_NONE) {
536 CeedScalar *e_vec_array;
537
538 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
539 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
540 }
541 }
542
543 // Q function
544 CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
545
546 // Restore input arrays
547 for (CeedInt i = 0; i < num_input_fields; i++) {
548 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, active_e_vec, false, impl));
549 }
550
551 // Output basis and restriction
552 for (CeedInt i = 0; i < num_output_fields; i++) {
553 bool is_active = false;
554 CeedInt field = impl->output_field_order[i];
555 CeedEvalMode eval_mode;
556 CeedVector l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field];
557
558 // Output vector
559 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec));
560 is_active = l_vec == CEED_VECTOR_ACTIVE;
561 if (is_active) {
562 l_vec = out_vec;
563 if (!e_vec) e_vec = active_e_vec;
564 }
565
566 // Basis action
567 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode));
568 switch (eval_mode) {
569 case CEED_EVAL_NONE:
570 break; // No action
571 case CEED_EVAL_INTERP:
572 case CEED_EVAL_GRAD:
573 case CEED_EVAL_DIV:
574 case CEED_EVAL_CURL: {
575 CeedBasis basis;
576
577 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis));
578 if (impl->apply_add_basis_out[field]) {
579 CeedCallBackend(CeedBasisApplyAdd(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec));
580 } else {
581 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec));
582 }
583 CeedCallBackend(CeedBasisDestroy(&basis));
584 break;
585 }
586 // LCOV_EXCL_START
587 case CEED_EVAL_WEIGHT: {
588 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
589 // LCOV_EXCL_STOP
590 }
591 }
592
593 // Restore evec
594 if (eval_mode == CEED_EVAL_NONE) {
595 CeedScalar *e_vec_array;
596
597 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[field], CEED_MEM_DEVICE, &e_vec_array));
598 CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
599 }
600
601 // Restrict
602 if (!impl->skip_rstr_out[field]) {
603 CeedElemRestriction elem_rstr;
604
605 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr));
606 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request));
607 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
608 }
609 if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
610 }
611
612 // Return work vector
613 CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec));
614 CeedCallBackend(CeedDestroy(&ceed));
615 CeedCallBackend(CeedQFunctionDestroy(&qf));
616 return CEED_ERROR_SUCCESS;
617 }
618
619 //------------------------------------------------------------------------------
620 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
621 //------------------------------------------------------------------------------
CeedOperatorSetupAtPoints_Cuda(CeedOperator op)622 static int CeedOperatorSetupAtPoints_Cuda(CeedOperator op) {
623 bool is_setup_done;
624 CeedInt max_num_points = -1, num_elem, num_input_fields, num_output_fields;
625 CeedQFunctionField *qf_input_fields, *qf_output_fields;
626 CeedQFunction qf;
627 CeedOperatorField *op_input_fields, *op_output_fields;
628 CeedOperator_Cuda *impl;
629
630 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
631 if (is_setup_done) return CEED_ERROR_SUCCESS;
632
633 CeedCallBackend(CeedOperatorGetData(op, &impl));
634 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
635 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
636 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
637 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
638 {
639 CeedElemRestriction rstr_points = NULL;
640
641 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
642 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
643 CeedCallBackend(CeedCalloc(num_elem, &impl->num_points));
644 for (CeedInt e = 0; e < num_elem; e++) {
645 CeedInt num_points_elem;
646
647 CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
648 impl->num_points[e] = num_points_elem;
649 }
650 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
651 }
652 impl->max_num_points = max_num_points;
653
654 // Allocate
655 CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in));
656 CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out));
657 CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in));
658 CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out));
659 CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out));
660 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order));
661 CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order));
662 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states));
663 CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in));
664 CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out));
665 impl->num_inputs = num_input_fields;
666 impl->num_outputs = num_output_fields;
667
668 // Set up infield and outfield e-vecs and q-vecs
669 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, true, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields,
670 max_num_points, num_elem));
671 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out,
672 impl->q_vecs_out, num_output_fields, max_num_points, num_elem));
673
674 // Reorder fields to allow reuse of buffers
675 impl->max_active_e_vec_len = 0;
676 {
677 bool is_ordered[CEED_FIELD_MAX];
678 CeedInt curr_index = 0;
679
680 for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
681 for (CeedInt i = 0; i < num_input_fields; i++) {
682 CeedSize e_vec_len_i;
683 CeedVector vec_i;
684 CeedElemRestriction rstr_i;
685
686 if (is_ordered[i]) continue;
687 is_ordered[i] = true;
688 impl->input_field_order[curr_index] = i;
689 curr_index++;
690 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
691 if (vec_i == CEED_VECTOR_NONE) {
692 // CEED_EVAL_WEIGHT
693 CeedCallBackend(CeedVectorDestroy(&vec_i));
694 continue;
695 };
696 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
697 CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
698 impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
699 for (CeedInt j = i + 1; j < num_input_fields; j++) {
700 CeedVector vec_j;
701 CeedElemRestriction rstr_j;
702
703 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
704 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
705 if (rstr_i == rstr_j && vec_i == vec_j) {
706 is_ordered[j] = true;
707 impl->input_field_order[curr_index] = j;
708 curr_index++;
709 }
710 CeedCallBackend(CeedVectorDestroy(&vec_j));
711 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
712 }
713 CeedCallBackend(CeedVectorDestroy(&vec_i));
714 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
715 }
716 }
717 {
718 bool is_ordered[CEED_FIELD_MAX];
719 CeedInt curr_index = 0;
720
721 for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false;
722 for (CeedInt i = 0; i < num_output_fields; i++) {
723 CeedSize e_vec_len_i;
724 CeedVector vec_i;
725 CeedElemRestriction rstr_i;
726
727 if (is_ordered[i]) continue;
728 is_ordered[i] = true;
729 impl->output_field_order[curr_index] = i;
730 curr_index++;
731 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i));
732 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i));
733 CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
734 impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
735 for (CeedInt j = i + 1; j < num_output_fields; j++) {
736 CeedVector vec_j;
737 CeedElemRestriction rstr_j;
738
739 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j));
740 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j));
741 if (rstr_i == rstr_j && vec_i == vec_j) {
742 is_ordered[j] = true;
743 impl->output_field_order[curr_index] = j;
744 curr_index++;
745 }
746 CeedCallBackend(CeedVectorDestroy(&vec_j));
747 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
748 }
749 CeedCallBackend(CeedVectorDestroy(&vec_i));
750 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
751 }
752 }
753 CeedCallBackend(CeedClearWorkVectors(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len));
754 {
755 // Create two work vectors for diagonal assembly
756 CeedVector temp_1, temp_2;
757
758 CeedCallBackend(CeedGetWorkVector(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len, &temp_1));
759 CeedCallBackend(CeedGetWorkVector(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len, &temp_2));
760 CeedCallBackend(CeedRestoreWorkVector(CeedOperatorReturnCeed(op), &temp_1));
761 CeedCallBackend(CeedRestoreWorkVector(CeedOperatorReturnCeed(op), &temp_2));
762 }
763 CeedCallBackend(CeedOperatorSetSetupDone(op));
764 CeedCallBackend(CeedQFunctionDestroy(&qf));
765 return CEED_ERROR_SUCCESS;
766 }
767
768 //------------------------------------------------------------------------------
769 // Input Basis Action AtPoints
770 //------------------------------------------------------------------------------
CeedOperatorInputBasisAtPoints_Cuda(CeedOperatorField op_input_field,CeedQFunctionField qf_input_field,CeedInt input_field,CeedVector in_vec,CeedVector active_e_vec,CeedInt num_elem,const CeedInt * num_points,const bool skip_active,const bool skip_passive,CeedOperator_Cuda * impl)771 static inline int CeedOperatorInputBasisAtPoints_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
772 CeedVector in_vec, CeedVector active_e_vec, CeedInt num_elem, const CeedInt *num_points,
773 const bool skip_active, const bool skip_passive, CeedOperator_Cuda *impl) {
774 bool is_active = false;
775 CeedEvalMode eval_mode;
776 CeedVector l_vec, e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field];
777
778 // Skip active input
779 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
780 is_active = l_vec == CEED_VECTOR_ACTIVE;
781 if (skip_active && is_active) return CEED_ERROR_SUCCESS;
782 if (skip_passive && !is_active) {
783 CeedCallBackend(CeedVectorDestroy(&l_vec));
784 return CEED_ERROR_SUCCESS;
785 }
786 if (is_active) {
787 l_vec = in_vec;
788 if (!e_vec) e_vec = active_e_vec;
789 }
790
791 // Basis action
792 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
793 switch (eval_mode) {
794 case CEED_EVAL_NONE: {
795 const CeedScalar *e_vec_array;
796
797 if (e_vec) {
798 CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, &e_vec_array));
799 } else {
800 CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, &e_vec_array));
801 }
802 CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
803 break;
804 }
805 case CEED_EVAL_INTERP:
806 case CEED_EVAL_GRAD:
807 case CEED_EVAL_DIV:
808 case CEED_EVAL_CURL: {
809 CeedBasis basis;
810
811 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis));
812 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, e_vec, q_vec));
813 CeedCallBackend(CeedBasisDestroy(&basis));
814 break;
815 }
816 case CEED_EVAL_WEIGHT:
817 break; // No action
818 }
819 if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
820 return CEED_ERROR_SUCCESS;
821 }
822
823 //------------------------------------------------------------------------------
824 // Apply and add to output AtPoints
825 //------------------------------------------------------------------------------
CeedOperatorApplyAddAtPoints_Cuda(CeedOperator op,CeedVector in_vec,CeedVector out_vec,CeedRequest * request)826 static int CeedOperatorApplyAddAtPoints_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
827 CeedInt max_num_points, *num_points, num_elem, num_input_fields, num_output_fields;
828 Ceed ceed;
829 CeedVector active_e_vec;
830 CeedQFunctionField *qf_input_fields, *qf_output_fields;
831 CeedQFunction qf;
832 CeedOperatorField *op_input_fields, *op_output_fields;
833 CeedOperator_Cuda *impl;
834
835 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
836 CeedCallBackend(CeedOperatorGetData(op, &impl));
837 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
838 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
839 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
840 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
841
842 // Setup
843 CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op));
844 num_points = impl->num_points;
845 max_num_points = impl->max_num_points;
846
847 // Work vector
848 CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec));
849
850 // Get point coordinates
851 {
852 CeedVector point_coords = NULL;
853 CeedElemRestriction rstr_points = NULL;
854
855 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords));
856 if (!impl->point_coords_elem) CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem));
857 {
858 uint64_t state;
859 CeedCallBackend(CeedVectorGetState(point_coords, &state));
860 if (impl->points_state != state) {
861 CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request));
862 }
863 }
864 CeedCallBackend(CeedVectorDestroy(&point_coords));
865 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
866 }
867
868 // Process inputs
869 for (CeedInt i = 0; i < num_input_fields; i++) {
870 CeedInt field = impl->input_field_order[i];
871
872 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, false, impl,
873 request));
874 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, num_elem,
875 num_points, false, false, impl));
876 }
877
878 // Output pointers, as necessary
879 for (CeedInt i = 0; i < num_output_fields; i++) {
880 CeedEvalMode eval_mode;
881
882 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
883 if (eval_mode == CEED_EVAL_NONE) {
884 CeedScalar *e_vec_array;
885
886 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
887 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
888 }
889 }
890
891 // Q function
892 CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out));
893
894 // Restore input arrays
895 for (CeedInt i = 0; i < num_input_fields; i++) {
896 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, active_e_vec, false, impl));
897 }
898
899 // Output basis and restriction
900 for (CeedInt i = 0; i < num_output_fields; i++) {
901 bool is_active = false;
902 CeedInt field = impl->output_field_order[i];
903 CeedEvalMode eval_mode;
904 CeedVector l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field];
905
906 // Output vector
907 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec));
908 is_active = l_vec == CEED_VECTOR_ACTIVE;
909 if (is_active) {
910 l_vec = out_vec;
911 if (!e_vec) e_vec = active_e_vec;
912 }
913
914 // Basis action
915 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode));
916 switch (eval_mode) {
917 case CEED_EVAL_NONE:
918 break; // No action
919 case CEED_EVAL_INTERP:
920 case CEED_EVAL_GRAD:
921 case CEED_EVAL_DIV:
922 case CEED_EVAL_CURL: {
923 CeedBasis basis;
924
925 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis));
926 if (impl->apply_add_basis_out[field]) {
927 CeedCallBackend(CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
928 } else {
929 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
930 }
931 CeedCallBackend(CeedBasisDestroy(&basis));
932 break;
933 }
934 // LCOV_EXCL_START
935 case CEED_EVAL_WEIGHT: {
936 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
937 // LCOV_EXCL_STOP
938 }
939 }
940
941 // Restore evec
942 if (eval_mode == CEED_EVAL_NONE) {
943 CeedScalar *e_vec_array;
944
945 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[field], CEED_MEM_DEVICE, &e_vec_array));
946 CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
947 }
948
949 // Restrict
950 if (!impl->skip_rstr_out[field]) {
951 CeedElemRestriction elem_rstr;
952
953 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr));
954 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request));
955 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
956 }
957 if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
958 }
959
960 // Restore work vector
961 CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec));
962 CeedCallBackend(CeedDestroy(&ceed));
963 CeedCallBackend(CeedQFunctionDestroy(&qf));
964 return CEED_ERROR_SUCCESS;
965 }
966
967 //------------------------------------------------------------------------------
968 // Linear QFunction Assembly Core
969 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op,bool build_objects,CeedVector * assembled,CeedElemRestriction * rstr,CeedRequest * request)970 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
971 CeedRequest *request) {
972 Ceed ceed, ceed_parent;
973 CeedInt num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
974 CeedScalar *assembled_array;
975 CeedVector *active_inputs;
976 CeedQFunctionField *qf_input_fields, *qf_output_fields;
977 CeedQFunction qf;
978 CeedOperatorField *op_input_fields, *op_output_fields;
979 CeedOperator_Cuda *impl;
980
981 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
982 CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
983 CeedCallBackend(CeedOperatorGetData(op, &impl));
984 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
985 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
986 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
987 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
988 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
989 active_inputs = impl->qf_active_in;
990 num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
991
992 // Setup
993 CeedCallBackend(CeedOperatorSetup_Cuda(op));
994
995 // Process inputs
996 for (CeedInt i = 0; i < num_input_fields; i++) {
997 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl, request));
998 CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, num_elem, true, impl));
999 }
1000
1001 // Count number of active input fields
1002 if (!num_active_in) {
1003 for (CeedInt i = 0; i < num_input_fields; i++) {
1004 CeedScalar *q_vec_array;
1005 CeedVector l_vec;
1006
1007 // Check if active input
1008 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec));
1009 if (l_vec == CEED_VECTOR_ACTIVE) {
1010 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
1011 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
1012 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
1013 CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs));
1014 for (CeedInt field = 0; field < size; field++) {
1015 CeedSize q_size = (CeedSize)Q * num_elem;
1016
1017 CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field]));
1018 CeedCallBackend(CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER,
1019 &q_vec_array[field * Q * num_elem]));
1020 }
1021 num_active_in += size;
1022 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
1023 }
1024 CeedCallBackend(CeedVectorDestroy(&l_vec));
1025 }
1026 impl->num_active_in = num_active_in;
1027 impl->qf_active_in = active_inputs;
1028 }
1029
1030 // Count number of active output fields
1031 if (!num_active_out) {
1032 for (CeedInt i = 0; i < num_output_fields; i++) {
1033 CeedVector l_vec;
1034
1035 // Check if active output
1036 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &l_vec));
1037 if (l_vec == CEED_VECTOR_ACTIVE) {
1038 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
1039 num_active_out += size;
1040 }
1041 CeedCallBackend(CeedVectorDestroy(&l_vec));
1042 }
1043 impl->num_active_out = num_active_out;
1044 }
1045
1046 // Check sizes
1047 CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
1048
1049 // Build objects if needed
1050 if (build_objects) {
1051 CeedSize l_size = (CeedSize)num_elem * Q * num_active_in * num_active_out;
1052 CeedInt strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
1053
1054 // Create output restriction
1055 CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
1056 (CeedSize)num_active_in * (CeedSize)num_active_out * (CeedSize)num_elem * (CeedSize)Q, strides,
1057 rstr));
1058 // Create assembled vector
1059 CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
1060 }
1061 CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
1062 CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
1063
1064 // Assemble QFunction
1065 for (CeedInt in = 0; in < num_active_in; in++) {
1066 // Set Inputs
1067 CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0));
1068 if (num_active_in > 1) {
1069 CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0));
1070 }
1071 // Set Outputs
1072 for (CeedInt out = 0; out < num_output_fields; out++) {
1073 CeedVector l_vec;
1074
1075 // Check if active output
1076 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec));
1077 if (l_vec == CEED_VECTOR_ACTIVE) {
1078 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
1079 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
1080 assembled_array += size * Q * num_elem; // Advance the pointer by the size of the output
1081 }
1082 CeedCallBackend(CeedVectorDestroy(&l_vec));
1083 }
1084 // Apply QFunction
1085 CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
1086 }
1087
1088 // Un-set output q-vecs to prevent accidental overwrite of Assembled
1089 for (CeedInt out = 0; out < num_output_fields; out++) {
1090 CeedVector l_vec;
1091
1092 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec));
1093 if (l_vec == CEED_VECTOR_ACTIVE) {
1094 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
1095 }
1096 CeedCallBackend(CeedVectorDestroy(&l_vec));
1097 }
1098
1099 // Restore input arrays
1100 for (CeedInt i = 0; i < num_input_fields; i++) {
1101 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl));
1102 }
1103
1104 // Restore output
1105 CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
1106 CeedCallBackend(CeedDestroy(&ceed));
1107 CeedCallBackend(CeedDestroy(&ceed_parent));
1108 CeedCallBackend(CeedQFunctionDestroy(&qf));
1109 return CEED_ERROR_SUCCESS;
1110 }
1111
1112 //------------------------------------------------------------------------------
1113 // Assemble Linear QFunction
1114 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op,CeedVector * assembled,CeedElemRestriction * rstr,CeedRequest * request)1115 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1116 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
1117 }
1118
1119 //------------------------------------------------------------------------------
1120 // Update Assembled Linear QFunction
1121 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op,CeedVector assembled,CeedElemRestriction rstr,CeedRequest * request)1122 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
1123 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
1124 }
1125
1126 //------------------------------------------------------------------------------
1127 // Assemble Diagonal Setup
1128 //------------------------------------------------------------------------------
CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op)1129 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) {
1130 Ceed ceed;
1131 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1132 CeedInt q_comp, num_nodes, num_qpts;
1133 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL;
1134 CeedBasis basis_in = NULL, basis_out = NULL;
1135 CeedQFunctionField *qf_fields;
1136 CeedQFunction qf;
1137 CeedOperatorField *op_fields;
1138 CeedOperator_Cuda *impl;
1139
1140 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1141 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1142 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
1143
1144 // Determine active input basis
1145 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1146 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1147 for (CeedInt i = 0; i < num_input_fields; i++) {
1148 CeedVector vec;
1149
1150 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1151 if (vec == CEED_VECTOR_ACTIVE) {
1152 CeedEvalMode eval_mode;
1153 CeedBasis basis;
1154
1155 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1156 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND,
1157 "Backend does not implement operator diagonal assembly with multiple active bases");
1158 if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1159 CeedCallBackend(CeedBasisDestroy(&basis));
1160 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1161 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1162 if (eval_mode != CEED_EVAL_WEIGHT) {
1163 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
1164 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1165 for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode;
1166 num_eval_modes_in += q_comp;
1167 }
1168 }
1169 CeedCallBackend(CeedVectorDestroy(&vec));
1170 }
1171
1172 // Determine active output basis
1173 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1174 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1175 for (CeedInt i = 0; i < num_output_fields; i++) {
1176 CeedVector vec;
1177
1178 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1179 if (vec == CEED_VECTOR_ACTIVE) {
1180 CeedBasis basis;
1181 CeedEvalMode eval_mode;
1182
1183 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1184 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1185 "Backend does not implement operator diagonal assembly with multiple active bases");
1186 if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1187 CeedCallBackend(CeedBasisDestroy(&basis));
1188 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1189 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1190 if (eval_mode != CEED_EVAL_WEIGHT) {
1191 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
1192 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1193 for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode;
1194 num_eval_modes_out += q_comp;
1195 }
1196 }
1197 CeedCallBackend(CeedVectorDestroy(&vec));
1198 }
1199
1200 // Operator data struct
1201 CeedCallBackend(CeedOperatorGetData(op, &impl));
1202 CeedCallBackend(CeedCalloc(1, &impl->diag));
1203 CeedOperatorDiag_Cuda *diag = impl->diag;
1204
1205 // Basis matrices
1206 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
1207 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
1208 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1209 const CeedInt interp_bytes = num_nodes * num_qpts * sizeof(CeedScalar);
1210 const CeedInt eval_modes_bytes = sizeof(CeedEvalMode);
1211 bool has_eval_none = false;
1212
1213 // CEED_EVAL_NONE
1214 for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1215 for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1216 if (has_eval_none) {
1217 CeedScalar *identity = NULL;
1218
1219 CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity));
1220 for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
1221 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes));
1222 CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice));
1223 CeedCallBackend(CeedFree(&identity));
1224 }
1225
1226 // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL
1227 for (CeedInt in = 0; in < 2; in++) {
1228 CeedFESpace fespace;
1229 CeedBasis basis = in ? basis_in : basis_out;
1230
1231 CeedCallBackend(CeedBasisGetFESpace(basis, &fespace));
1232 switch (fespace) {
1233 case CEED_FE_SPACE_H1: {
1234 CeedInt q_comp_interp, q_comp_grad;
1235 const CeedScalar *interp, *grad;
1236 CeedScalar *d_interp, *d_grad;
1237
1238 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1239 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
1240
1241 CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1242 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1243 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1244 CeedCallBackend(CeedBasisGetGrad(basis, &grad));
1245 CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad));
1246 CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice));
1247 if (in) {
1248 diag->d_interp_in = d_interp;
1249 diag->d_grad_in = d_grad;
1250 } else {
1251 diag->d_interp_out = d_interp;
1252 diag->d_grad_out = d_grad;
1253 }
1254 } break;
1255 case CEED_FE_SPACE_HDIV: {
1256 CeedInt q_comp_interp, q_comp_div;
1257 const CeedScalar *interp, *div;
1258 CeedScalar *d_interp, *d_div;
1259
1260 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1261 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
1262
1263 CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1264 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1265 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1266 CeedCallBackend(CeedBasisGetDiv(basis, &div));
1267 CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div));
1268 CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice));
1269 if (in) {
1270 diag->d_interp_in = d_interp;
1271 diag->d_div_in = d_div;
1272 } else {
1273 diag->d_interp_out = d_interp;
1274 diag->d_div_out = d_div;
1275 }
1276 } break;
1277 case CEED_FE_SPACE_HCURL: {
1278 CeedInt q_comp_interp, q_comp_curl;
1279 const CeedScalar *interp, *curl;
1280 CeedScalar *d_interp, *d_curl;
1281
1282 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1283 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
1284
1285 CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1286 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1287 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1288 CeedCallBackend(CeedBasisGetCurl(basis, &curl));
1289 CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl));
1290 CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice));
1291 if (in) {
1292 diag->d_interp_in = d_interp;
1293 diag->d_curl_in = d_curl;
1294 } else {
1295 diag->d_interp_out = d_interp;
1296 diag->d_curl_out = d_curl;
1297 }
1298 } break;
1299 }
1300 }
1301
1302 // Arrays of eval_modes
1303 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes));
1304 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice));
1305 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes));
1306 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice));
1307 CeedCallBackend(CeedFree(&eval_modes_in));
1308 CeedCallBackend(CeedFree(&eval_modes_out));
1309 CeedCallBackend(CeedDestroy(&ceed));
1310 CeedCallBackend(CeedBasisDestroy(&basis_in));
1311 CeedCallBackend(CeedBasisDestroy(&basis_out));
1312 CeedCallBackend(CeedQFunctionDestroy(&qf));
1313 return CEED_ERROR_SUCCESS;
1314 }
1315
1316 //------------------------------------------------------------------------------
1317 // Assemble Diagonal Setup (Compilation)
1318 //------------------------------------------------------------------------------
CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op,CeedInt use_ceedsize_idx,const bool is_point_block)1319 static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) {
1320 Ceed ceed;
1321 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1322 CeedInt num_comp, q_comp, num_nodes, num_qpts;
1323 CeedBasis basis_in = NULL, basis_out = NULL;
1324 CeedQFunctionField *qf_fields;
1325 CeedQFunction qf;
1326 CeedOperatorField *op_fields;
1327 CeedOperator_Cuda *impl;
1328
1329 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1330 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1331 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
1332
1333 // Determine active input basis
1334 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1335 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1336 for (CeedInt i = 0; i < num_input_fields; i++) {
1337 CeedVector vec;
1338
1339 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1340 if (vec == CEED_VECTOR_ACTIVE) {
1341 CeedEvalMode eval_mode;
1342 CeedBasis basis;
1343
1344 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1345 if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1346 CeedCallBackend(CeedBasisDestroy(&basis));
1347 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1348 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1349 if (eval_mode != CEED_EVAL_WEIGHT) {
1350 num_eval_modes_in += q_comp;
1351 }
1352 }
1353 CeedCallBackend(CeedVectorDestroy(&vec));
1354 }
1355
1356 // Determine active output basis
1357 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1358 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1359 for (CeedInt i = 0; i < num_output_fields; i++) {
1360 CeedVector vec;
1361
1362 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1363 if (vec == CEED_VECTOR_ACTIVE) {
1364 CeedEvalMode eval_mode;
1365 CeedBasis basis;
1366
1367 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1368 if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1369 CeedCallBackend(CeedBasisDestroy(&basis));
1370 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1371 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1372 if (eval_mode != CEED_EVAL_WEIGHT) {
1373 num_eval_modes_out += q_comp;
1374 }
1375 }
1376 CeedCallBackend(CeedVectorDestroy(&vec));
1377 }
1378
1379 // Operator data struct
1380 CeedCallBackend(CeedOperatorGetData(op, &impl));
1381 CeedOperatorDiag_Cuda *diag = impl->diag;
1382
1383 // Assemble kernel
1384 const char diagonal_kernel_source[] = "// Diagonal assembly source\n#include <ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h>\n";
1385 CUmodule *module = is_point_block ? &diag->module_point_block : &diag->module;
1386 CeedInt elems_per_block = 1;
1387
1388 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
1389 CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
1390 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
1391 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1392 CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1393 num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE",
1394 use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block));
1395 CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal));
1396 CeedCallBackend(CeedDestroy(&ceed));
1397 CeedCallBackend(CeedBasisDestroy(&basis_in));
1398 CeedCallBackend(CeedBasisDestroy(&basis_out));
1399 CeedCallBackend(CeedQFunctionDestroy(&qf));
1400 return CEED_ERROR_SUCCESS;
1401 }
1402
1403 //------------------------------------------------------------------------------
1404 // Assemble Diagonal Core
1405 //------------------------------------------------------------------------------
CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op,CeedVector assembled,CeedRequest * request,const bool is_point_block)1406 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
1407 Ceed ceed;
1408 CeedInt num_elem, num_nodes;
1409 CeedScalar *elem_diag_array;
1410 const CeedScalar *assembled_qf_array;
1411 CeedVector assembled_qf = NULL, elem_diag;
1412 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr;
1413 CeedOperator_Cuda *impl;
1414
1415 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1416 CeedCallBackend(CeedOperatorGetData(op, &impl));
1417
1418 // Assemble QFunction
1419 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request));
1420 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1421 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1422
1423 // Setup
1424 if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op));
1425 CeedOperatorDiag_Cuda *diag = impl->diag;
1426
1427 assert(diag != NULL);
1428
1429 // Assemble kernel if needed
1430 if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) {
1431 CeedSize assembled_length, assembled_qf_length;
1432 CeedInt use_ceedsize_idx = 0;
1433 CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
1434 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1435 if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1436
1437 CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block));
1438 }
1439
1440 // Restriction and diagonal vector
1441 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1442 CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND,
1443 "Cannot assemble operator diagonal with different input and output active element restrictions");
1444 if (!is_point_block && !diag->diag_rstr) {
1445 CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr));
1446 CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag));
1447 } else if (is_point_block && !diag->point_block_diag_rstr) {
1448 CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr));
1449 CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag));
1450 }
1451 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1452 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1453 diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
1454 elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
1455 CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
1456
1457 // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers
1458 CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes));
1459 if (num_nodes > 0) {
1460 // Assemble element operator diagonals
1461 CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
1462 CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
1463
1464 // Compute the diagonal of B^T D B
1465 CeedInt elems_per_block = 1;
1466 CeedInt grid = CeedDivUpInt(num_elem, elems_per_block);
1467 void *args[] = {(void *)&num_elem, &diag->d_identity, &diag->d_interp_in, &diag->d_grad_in, &diag->d_div_in,
1468 &diag->d_curl_in, &diag->d_interp_out, &diag->d_grad_out, &diag->d_div_out, &diag->d_curl_out,
1469 &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array};
1470
1471 if (is_point_block) {
1472 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args));
1473 } else {
1474 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args));
1475 }
1476
1477 // Restore arrays
1478 CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
1479 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1480 }
1481
1482 // Assemble local operator diagonal
1483 CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
1484
1485 // Cleanup
1486 CeedCallBackend(CeedDestroy(&ceed));
1487 CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1488 return CEED_ERROR_SUCCESS;
1489 }
1490
1491 //------------------------------------------------------------------------------
1492 // Assemble Linear Diagonal
1493 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op,CeedVector assembled,CeedRequest * request)1494 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1495 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
1496 return CEED_ERROR_SUCCESS;
1497 }
1498
1499 //------------------------------------------------------------------------------
1500 // Assemble Linear Point Block Diagonal
1501 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op,CeedVector assembled,CeedRequest * request)1502 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1503 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
1504 return CEED_ERROR_SUCCESS;
1505 }
1506
1507 //------------------------------------------------------------------------------
1508 // Single Operator Assembly Setup
1509 //------------------------------------------------------------------------------
CeedOperatorAssembleSingleSetup_Cuda(CeedOperator op,CeedInt use_ceedsize_idx)1510 static int CeedOperatorAssembleSingleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) {
1511 Ceed ceed;
1512 Ceed_Cuda *cuda_data;
1513 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1514 CeedInt elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp;
1515 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL;
1516 CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
1517 CeedBasis basis_in = NULL, basis_out = NULL;
1518 CeedQFunctionField *qf_fields;
1519 CeedQFunction qf;
1520 CeedOperatorField *input_fields, *output_fields;
1521 CeedOperator_Cuda *impl;
1522
1523 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1524 CeedCallBackend(CeedOperatorGetData(op, &impl));
1525
1526 // Get intput and output fields
1527 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
1528
1529 // Determine active input basis eval mode
1530 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1531 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1532 for (CeedInt i = 0; i < num_input_fields; i++) {
1533 CeedVector vec;
1534
1535 CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
1536 if (vec == CEED_VECTOR_ACTIVE) {
1537 CeedEvalMode eval_mode;
1538 CeedElemRestriction elem_rstr;
1539 CeedBasis basis;
1540
1541 CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis));
1542 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases");
1543 if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1544 CeedCallBackend(CeedBasisDestroy(&basis));
1545 CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &elem_rstr));
1546 if (!rstr_in) CeedCallBackend(CeedElemRestrictionReferenceCopy(elem_rstr, &rstr_in));
1547 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1548 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1549 if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in;
1550 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in));
1551 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1552 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1553 if (eval_mode != CEED_EVAL_WEIGHT) {
1554 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1555 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1556 for (CeedInt d = 0; d < q_comp; d++) {
1557 eval_modes_in[num_eval_modes_in + d] = eval_mode;
1558 }
1559 num_eval_modes_in += q_comp;
1560 }
1561 }
1562 CeedCallBackend(CeedVectorDestroy(&vec));
1563 }
1564
1565 // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
1566 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1567 for (CeedInt i = 0; i < num_output_fields; i++) {
1568 CeedVector vec;
1569
1570 CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1571 if (vec == CEED_VECTOR_ACTIVE) {
1572 CeedEvalMode eval_mode;
1573 CeedElemRestriction elem_rstr;
1574 CeedBasis basis;
1575
1576 CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis));
1577 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1578 "Backend does not implement operator assembly with multiple active bases");
1579 if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1580 CeedCallBackend(CeedBasisDestroy(&basis));
1581 CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &elem_rstr));
1582 if (!rstr_out) CeedCallBackend(CeedElemRestrictionReferenceCopy(elem_rstr, &rstr_out));
1583 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1584 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1585 if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out;
1586 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out));
1587 CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED,
1588 "Active input and output bases must have the same number of quadrature points");
1589 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1590 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1591 if (eval_mode != CEED_EVAL_WEIGHT) {
1592 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1593 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1594 for (CeedInt d = 0; d < q_comp; d++) {
1595 eval_modes_out[num_eval_modes_out + d] = eval_mode;
1596 }
1597 num_eval_modes_out += q_comp;
1598 }
1599 }
1600 CeedCallBackend(CeedVectorDestroy(&vec));
1601 }
1602 CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1603
1604 CeedCallBackend(CeedCalloc(1, &impl->asmb));
1605 CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1606 asmb->elems_per_block = 1;
1607 asmb->block_size_x = elem_size_in;
1608 asmb->block_size_y = elem_size_out;
1609
1610 CeedCallBackend(CeedGetData(ceed, &cuda_data));
1611 bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock;
1612
1613 if (fallback) {
1614 // Use fallback kernel with 1D threadblock
1615 asmb->block_size_y = 1;
1616 }
1617
1618 // Compile kernels
1619 const char assembly_kernel_source[] = "// Full assembly source\n#include <ceed/jit-source/cuda/cuda-ref-operator-assemble.h>\n";
1620
1621 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in));
1622 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out));
1623 CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1624 num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in,
1625 "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE",
1626 asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y,
1627 "USE_CEEDSIZE", use_ceedsize_idx));
1628 CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble));
1629
1630 // Load into B_in, in order that they will be used in eval_modes_in
1631 {
1632 const CeedInt in_bytes = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar);
1633 CeedInt d_in = 0;
1634 CeedEvalMode eval_modes_in_prev = CEED_EVAL_NONE;
1635 bool has_eval_none = false;
1636 CeedScalar *identity = NULL;
1637
1638 for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1639 has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1640 }
1641 if (has_eval_none) {
1642 CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity));
1643 for (CeedInt i = 0; i < (elem_size_in < num_qpts_in ? elem_size_in : num_qpts_in); i++) identity[i * elem_size_in + i] = 1.0;
1644 }
1645
1646 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes));
1647 for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1648 const CeedScalar *h_B_in;
1649
1650 CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in));
1651 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp));
1652 if (q_comp > 1) {
1653 if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0;
1654 else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in];
1655 }
1656 eval_modes_in_prev = eval_modes_in[i];
1657
1658 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[i * elem_size_in * num_qpts_in], h_B_in, elem_size_in * num_qpts_in * sizeof(CeedScalar),
1659 cudaMemcpyHostToDevice));
1660 }
1661 CeedCallBackend(CeedFree(&identity));
1662 }
1663 CeedCallBackend(CeedFree(&eval_modes_in));
1664
1665 // Load into B_out, in order that they will be used in eval_modes_out
1666 {
1667 const CeedInt out_bytes = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar);
1668 CeedInt d_out = 0;
1669 CeedEvalMode eval_modes_out_prev = CEED_EVAL_NONE;
1670 bool has_eval_none = false;
1671 CeedScalar *identity = NULL;
1672
1673 for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1674 has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1675 }
1676 if (has_eval_none) {
1677 CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity));
1678 for (CeedInt i = 0; i < (elem_size_out < num_qpts_out ? elem_size_out : num_qpts_out); i++) identity[i * elem_size_out + i] = 1.0;
1679 }
1680
1681 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes));
1682 for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1683 const CeedScalar *h_B_out;
1684
1685 CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out));
1686 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp));
1687 if (q_comp > 1) {
1688 if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0;
1689 else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out];
1690 }
1691 eval_modes_out_prev = eval_modes_out[i];
1692
1693 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[i * elem_size_out * num_qpts_out], h_B_out, elem_size_out * num_qpts_out * sizeof(CeedScalar),
1694 cudaMemcpyHostToDevice));
1695 }
1696 CeedCallBackend(CeedFree(&identity));
1697 }
1698 CeedCallBackend(CeedFree(&eval_modes_out));
1699 CeedCallBackend(CeedDestroy(&ceed));
1700 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1701 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1702 CeedCallBackend(CeedBasisDestroy(&basis_in));
1703 CeedCallBackend(CeedBasisDestroy(&basis_out));
1704 CeedCallBackend(CeedQFunctionDestroy(&qf));
1705 return CEED_ERROR_SUCCESS;
1706 }
1707
1708 //------------------------------------------------------------------------------
1709 // Assemble matrix data for COO matrix of assembled operator.
1710 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1711 //
1712 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator
1713 // (could have multiple basis eval modes).
1714 // TODO: allow multiple active input restrictions/basis objects
1715 //------------------------------------------------------------------------------
CeedOperatorAssembleSingle_Cuda(CeedOperator op,CeedInt offset,CeedVector values)1716 static int CeedOperatorAssembleSingle_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1717 Ceed ceed;
1718 CeedSize values_length = 0, assembled_qf_length = 0;
1719 CeedInt use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out;
1720 CeedScalar *values_array;
1721 const CeedScalar *assembled_qf_array;
1722 CeedVector assembled_qf = NULL;
1723 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out;
1724 CeedRestrictionType rstr_type_in, rstr_type_out;
1725 const bool *orients_in = NULL, *orients_out = NULL;
1726 const CeedInt8 *curl_orients_in = NULL, *curl_orients_out = NULL;
1727 CeedOperator_Cuda *impl;
1728
1729 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1730 CeedCallBackend(CeedOperatorGetData(op, &impl));
1731
1732 // Assemble QFunction
1733 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE));
1734 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1735 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1736
1737 CeedCallBackend(CeedVectorGetLength(values, &values_length));
1738 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1739 if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1740
1741 // Setup
1742 if (!impl->asmb) CeedCallBackend(CeedOperatorAssembleSingleSetup_Cuda(op, use_ceedsize_idx));
1743 CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1744
1745 assert(asmb != NULL);
1746
1747 // Assemble element operator
1748 CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1749 values_array += offset;
1750
1751 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1752 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in));
1753 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1754
1755 CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in));
1756 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1757 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in));
1758 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1759 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in));
1760 }
1761
1762 if (rstr_in != rstr_out) {
1763 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out));
1764 CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED,
1765 "Active input and output operator restrictions must have the same number of elements");
1766 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1767
1768 CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out));
1769 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1770 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out));
1771 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1772 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out));
1773 }
1774 } else {
1775 elem_size_out = elem_size_in;
1776 orients_out = orients_in;
1777 curl_orients_out = curl_orients_in;
1778 }
1779
1780 // Compute B^T D B
1781 CeedInt shared_mem =
1782 ((curl_orients_in || curl_orients_out ? elem_size_in * elem_size_out : 0) + (curl_orients_in ? elem_size_in * asmb->block_size_y : 0)) *
1783 sizeof(CeedScalar);
1784 CeedInt grid = CeedDivUpInt(num_elem_in, asmb->elems_per_block);
1785 void *args[] = {(void *)&num_elem_in, &asmb->d_B_in, &asmb->d_B_out, &orients_in, &curl_orients_in,
1786 &orients_out, &curl_orients_out, &assembled_qf_array, &values_array};
1787
1788 CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, NULL, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block,
1789 shared_mem, args));
1790
1791 // Restore arrays
1792 CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1793 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1794
1795 // Cleanup
1796 CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1797 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1798 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in));
1799 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1800 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in));
1801 }
1802 if (rstr_in != rstr_out) {
1803 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1804 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out));
1805 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1806 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out));
1807 }
1808 }
1809 CeedCallBackend(CeedDestroy(&ceed));
1810 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1811 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1812 return CEED_ERROR_SUCCESS;
1813 }
1814
1815 //------------------------------------------------------------------------------
1816 // Assemble Linear QFunction AtPoints
1817 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleQFunctionAtPoints_Cuda(CeedOperator op,CeedVector * assembled,CeedElemRestriction * rstr,CeedRequest * request)1818 static int CeedOperatorLinearAssembleQFunctionAtPoints_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1819 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "Backend does not implement CeedOperatorLinearAssembleQFunction");
1820 }
1821
1822 //------------------------------------------------------------------------------
1823 // Assemble Linear Diagonal AtPoints
1824 //------------------------------------------------------------------------------
CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda(CeedOperator op,CeedVector assembled,CeedRequest * request)1825 static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1826 CeedInt max_num_points, *num_points, num_elem, num_input_fields, num_output_fields;
1827 Ceed ceed;
1828 CeedVector active_e_vec_in, active_e_vec_out;
1829 CeedQFunctionField *qf_input_fields, *qf_output_fields;
1830 CeedQFunction qf;
1831 CeedOperatorField *op_input_fields, *op_output_fields;
1832 CeedOperator_Cuda *impl;
1833
1834 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1835 CeedCallBackend(CeedOperatorGetData(op, &impl));
1836 CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1837 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
1838 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1839 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1840
1841 // Setup
1842 CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op));
1843 num_points = impl->num_points;
1844 max_num_points = impl->max_num_points;
1845
1846 // Work vector
1847 CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec_in));
1848 CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec_out));
1849 {
1850 CeedSize length_in, length_out;
1851
1852 CeedCallBackend(CeedVectorGetLength(active_e_vec_in, &length_in));
1853 CeedCallBackend(CeedVectorGetLength(active_e_vec_out, &length_out));
1854 // Need input e_vec to be longer
1855 if (length_in < length_out) {
1856 CeedVector temp = active_e_vec_in;
1857
1858 active_e_vec_in = active_e_vec_out;
1859 active_e_vec_out = temp;
1860 }
1861 }
1862
1863 // Get point coordinates
1864 {
1865 CeedVector point_coords = NULL;
1866 CeedElemRestriction rstr_points = NULL;
1867
1868 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords));
1869 if (!impl->point_coords_elem) CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem));
1870 {
1871 uint64_t state;
1872 CeedCallBackend(CeedVectorGetState(point_coords, &state));
1873 if (impl->points_state != state) {
1874 CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request));
1875 }
1876 }
1877 CeedCallBackend(CeedVectorDestroy(&point_coords));
1878 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1879 }
1880
1881 // Process inputs
1882 for (CeedInt i = 0; i < num_input_fields; i++) {
1883 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl, request));
1884 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, num_elem, num_points, true, false,
1885 impl));
1886 }
1887
1888 // Output pointers, as necessary
1889 for (CeedInt i = 0; i < num_output_fields; i++) {
1890 CeedEvalMode eval_mode;
1891
1892 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1893 if (eval_mode == CEED_EVAL_NONE) {
1894 CeedScalar *e_vec_array;
1895
1896 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
1897 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
1898 }
1899 }
1900
1901 // Loop over active fields
1902 for (CeedInt i = 0; i < num_input_fields; i++) {
1903 bool is_active = false, is_active_at_points = true;
1904 CeedInt elem_size = 1, num_comp_active = 1, e_vec_size = 0, field_in = impl->input_field_order[i];
1905 CeedRestrictionType rstr_type;
1906 CeedVector l_vec;
1907 CeedElemRestriction elem_rstr;
1908
1909 // -- Skip non-active input
1910 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[field_in], &l_vec));
1911 is_active = l_vec == CEED_VECTOR_ACTIVE;
1912 CeedCallBackend(CeedVectorDestroy(&l_vec));
1913 if (!is_active || impl->skip_rstr_in[field_in]) continue;
1914
1915 // -- Get active restriction type
1916 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[field_in], &elem_rstr));
1917 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
1918 is_active_at_points = rstr_type == CEED_RESTRICTION_POINTS;
1919 if (!is_active_at_points) CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1920 else elem_size = max_num_points;
1921 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp_active));
1922 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1923
1924 e_vec_size = elem_size * num_comp_active;
1925 CeedCallBackend(CeedVectorSetValue(active_e_vec_in, 0.0));
1926 for (CeedInt s = 0; s < e_vec_size; s++) {
1927 CeedVector q_vec = impl->q_vecs_in[field_in];
1928
1929 // Update unit vector
1930 {
1931 // Note: E-vec strides are node * (1) + comp * (elem_size * num_elem) + elem * (elem_size)
1932 CeedInt node = (s - 1) % elem_size, comp = (s - 1) / elem_size;
1933 CeedSize start = node * 1 + comp * (elem_size * num_elem);
1934 CeedSize stop = (comp + 1) * (elem_size * num_elem);
1935
1936 if (s != 0) CeedCallBackend(CeedVectorSetValueStrided(active_e_vec_in, start, stop, elem_size, 0.0));
1937
1938 node = s % elem_size, comp = s / elem_size;
1939 start = node * 1 + comp * (elem_size * num_elem);
1940 stop = (comp + 1) * (elem_size * num_elem);
1941 CeedCallBackend(CeedVectorSetValueStrided(active_e_vec_in, start, stop, elem_size, 1.0));
1942 }
1943
1944 // Basis action
1945 for (CeedInt j = 0; j < num_input_fields; j++) {
1946 CeedInt field = impl->input_field_order[j];
1947
1948 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[field], qf_input_fields[field], field, NULL, active_e_vec_in, num_elem,
1949 num_points, false, true, impl));
1950 }
1951
1952 // Q function
1953 CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out));
1954
1955 // Output basis apply if needed
1956 for (CeedInt j = 0; j < num_output_fields; j++) {
1957 bool is_active = false;
1958 CeedInt elem_size = 0;
1959 CeedInt field_out = impl->output_field_order[j];
1960 CeedRestrictionType rstr_type;
1961 CeedEvalMode eval_mode;
1962 CeedVector l_vec, e_vec = impl->e_vecs_out[field_out], q_vec = impl->q_vecs_out[field_out];
1963 CeedElemRestriction elem_rstr;
1964
1965 // ---- Skip non-active output
1966 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field_out], &l_vec));
1967 is_active = l_vec == CEED_VECTOR_ACTIVE;
1968 CeedCallBackend(CeedVectorDestroy(&l_vec));
1969 if (!is_active) continue;
1970 if (!e_vec) e_vec = active_e_vec_out;
1971
1972 // ---- Check if elem size matches
1973 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field_out], &elem_rstr));
1974 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
1975 if (is_active_at_points && rstr_type != CEED_RESTRICTION_POINTS) continue;
1976 if (rstr_type == CEED_RESTRICTION_POINTS) {
1977 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(elem_rstr, &elem_size));
1978 } else {
1979 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1980 }
1981 {
1982 CeedInt num_comp = 0;
1983
1984 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1985 if (e_vec_size != num_comp * elem_size) continue;
1986 }
1987
1988 // Basis action
1989 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field_out], &eval_mode));
1990 switch (eval_mode) {
1991 case CEED_EVAL_NONE: {
1992 CeedScalar *e_vec_array;
1993
1994 CeedCallBackend(CeedVectorTakeArray(q_vec, CEED_MEM_DEVICE, &e_vec_array));
1995 CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
1996 break;
1997 }
1998 case CEED_EVAL_INTERP:
1999 case CEED_EVAL_GRAD:
2000 case CEED_EVAL_DIV:
2001 case CEED_EVAL_CURL: {
2002 CeedBasis basis;
2003
2004 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field_out], &basis));
2005 if (impl->apply_add_basis_out[field_out]) {
2006 CeedCallBackend(CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec,
2007 e_vec));
2008 } else {
2009 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
2010 }
2011 CeedCallBackend(CeedBasisDestroy(&basis));
2012 break;
2013 }
2014 // LCOV_EXCL_START
2015 case CEED_EVAL_WEIGHT: {
2016 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
2017 // LCOV_EXCL_STOP
2018 }
2019 }
2020
2021 // Mask output e-vec
2022 if (impl->skip_rstr_out[field_out]) {
2023 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2024 continue;
2025 }
2026 CeedCallBackend(CeedVectorPointwiseMult(e_vec, active_e_vec_in, e_vec));
2027
2028 // Restrict
2029 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, assembled, request));
2030 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2031
2032 // Reset q_vec for
2033 if (eval_mode == CEED_EVAL_NONE) {
2034 CeedScalar *e_vec_array;
2035
2036 CeedCallBackend(CeedVectorGetArrayWrite(e_vec, CEED_MEM_DEVICE, &e_vec_array));
2037 CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
2038 }
2039 }
2040
2041 // Reset vec
2042 if (s == e_vec_size - 1 && i != num_input_fields - 1) CeedCallBackend(CeedVectorSetValue(q_vec, 0.0));
2043 }
2044 }
2045
2046 // Restore CEED_EVAL_NONE
2047 for (CeedInt i = 0; i < num_output_fields; i++) {
2048 CeedEvalMode eval_mode;
2049
2050 // Get eval_mode
2051 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2052
2053 // Restore evec
2054 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2055 if (eval_mode == CEED_EVAL_NONE) {
2056 CeedScalar *e_vec_array;
2057
2058 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
2059 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs_out[i], &e_vec_array));
2060 }
2061 }
2062
2063 // Restore input arrays
2064 for (CeedInt i = 0; i < num_input_fields; i++) {
2065 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl));
2066 }
2067
2068 // Restore work vector
2069 CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec_in));
2070 CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec_out));
2071 CeedCallBackend(CeedDestroy(&ceed));
2072 CeedCallBackend(CeedQFunctionDestroy(&qf));
2073 return CEED_ERROR_SUCCESS;
2074 }
2075
2076 //------------------------------------------------------------------------------
2077 // Create operator
2078 //------------------------------------------------------------------------------
CeedOperatorCreate_Cuda(CeedOperator op)2079 int CeedOperatorCreate_Cuda(CeedOperator op) {
2080 Ceed ceed;
2081 CeedOperator_Cuda *impl;
2082
2083 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
2084 CeedCallBackend(CeedCalloc(1, &impl));
2085 CeedCallBackend(CeedOperatorSetData(op, impl));
2086
2087 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
2088 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
2089 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
2090 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal",
2091 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
2092 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedOperatorAssembleSingle_Cuda));
2093 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
2094 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
2095 CeedCallBackend(CeedDestroy(&ceed));
2096 return CEED_ERROR_SUCCESS;
2097 }
2098
2099 //------------------------------------------------------------------------------
2100 // Create operator AtPoints
2101 //------------------------------------------------------------------------------
CeedOperatorCreateAtPoints_Cuda(CeedOperator op)2102 int CeedOperatorCreateAtPoints_Cuda(CeedOperator op) {
2103 Ceed ceed;
2104 CeedOperator_Cuda *impl;
2105
2106 CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
2107 CeedCallBackend(CeedCalloc(1, &impl));
2108 CeedCallBackend(CeedOperatorSetData(op, impl));
2109
2110 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunctionAtPoints_Cuda));
2111 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda));
2112 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAddAtPoints_Cuda));
2113 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
2114 CeedCallBackend(CeedDestroy(&ceed));
2115 return CEED_ERROR_SUCCESS;
2116 }
2117
2118 //------------------------------------------------------------------------------
2119