xref: /libCEED/interface/ceed-preconditioning.c (revision bcd92680c8b93d5b8861c3c87761f161cf9f50ef)
1 // Copyright (c) 2017-2024, 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-impl.h>
9 #include <ceed.h>
10 #include <ceed/backend.h>
11 #include <assert.h>
12 #include <math.h>
13 #include <stdbool.h>
14 #include <stdio.h>
15 #include <string.h>
16 
17 /// @file
18 /// Implementation of CeedOperator preconditioning interfaces
19 
20 /// ----------------------------------------------------------------------------
21 /// CeedOperator Library Internal Preconditioning Functions
22 /// ----------------------------------------------------------------------------
23 /// @addtogroup CeedOperatorDeveloper
24 /// @{
25 
26 /**
27   @brief Duplicate a `CeedQFunction` with a reference `Ceed` to fallback for advanced `CeedOperator` functionality
28 
29   @param[in]  fallback_ceed `Ceed` on which to create fallback `CeedQFunction`
30   @param[in]  qf            `CeedQFunction` to create fallback for
31   @param[out] qf_fallback   Fallback `CeedQFunction`
32 
33   @return An error code: 0 - success, otherwise - failure
34 
35   @ref Developer
36 **/
37 static int CeedQFunctionCreateFallback(Ceed fallback_ceed, CeedQFunction qf, CeedQFunction *qf_fallback) {
38   char               *source_path_with_name = NULL;
39   CeedInt             num_input_fields, num_output_fields;
40   CeedQFunctionField *input_fields, *output_fields;
41 
42   // Check if NULL qf passed in
43   if (!qf) return CEED_ERROR_SUCCESS;
44 
45   CeedDebug256(CeedQFunctionReturnCeed(qf), 1, "---------- CeedOperator Fallback ----------\n");
46   CeedDebug(CeedQFunctionReturnCeed(qf), "Creating fallback CeedQFunction\n");
47 
48   if (qf->source_path) {
49     size_t path_len = strlen(qf->source_path), name_len = strlen(qf->kernel_name);
50 
51     CeedCall(CeedCalloc(path_len + name_len + 2, &source_path_with_name));
52     memcpy(source_path_with_name, qf->source_path, path_len);
53     memcpy(&source_path_with_name[path_len], ":", 1);
54     memcpy(&source_path_with_name[path_len + 1], qf->kernel_name, name_len);
55   } else if (qf->user_source) {
56     CeedCall(CeedStringAllocCopy(qf->user_source, &source_path_with_name));
57   } else {
58     CeedCall(CeedCalloc(1, &source_path_with_name));
59   }
60 
61   {
62     CeedInt           vec_length;
63     CeedQFunctionUser f;
64 
65     CeedCall(CeedQFunctionGetVectorLength(qf, &vec_length));
66     CeedCall(CeedQFunctionGetUserFunction(qf, &f));
67     CeedCall(CeedQFunctionCreateInterior(fallback_ceed, vec_length, f, source_path_with_name, qf_fallback));
68   }
69   {
70     CeedQFunctionContext ctx;
71 
72     CeedCall(CeedQFunctionGetContext(qf, &ctx));
73     CeedCall(CeedQFunctionSetContext(*qf_fallback, ctx));
74     CeedCall(CeedQFunctionContextDestroy(&ctx));
75   }
76   CeedCall(CeedQFunctionGetFields(qf, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
77   for (CeedInt i = 0; i < num_input_fields; i++) {
78     const char  *field_name;
79     CeedInt      size;
80     CeedEvalMode eval_mode;
81 
82     CeedCall(CeedQFunctionFieldGetData(input_fields[i], &field_name, &size, &eval_mode));
83     CeedCall(CeedQFunctionAddInput(*qf_fallback, field_name, size, eval_mode));
84   }
85   for (CeedInt i = 0; i < num_output_fields; i++) {
86     const char  *field_name;
87     CeedInt      size;
88     CeedEvalMode eval_mode;
89 
90     CeedCall(CeedQFunctionFieldGetData(output_fields[i], &field_name, &size, &eval_mode));
91     CeedCall(CeedQFunctionAddOutput(*qf_fallback, field_name, size, eval_mode));
92   }
93   CeedCall(CeedFree(&source_path_with_name));
94   return CEED_ERROR_SUCCESS;
95 }
96 
97 /**
98   @brief Duplicate a `CeedOperator` with a reference `Ceed` to fallback for advanced `CeedOperator` functionality
99 
100   @param[in,out] op `CeedOperator` to create fallback for
101 
102   @return An error code: 0 - success, otherwise - failure
103 
104   @ref Developer
105 **/
106 static int CeedOperatorCreateFallback(CeedOperator op) {
107   bool         is_composite;
108   Ceed         ceed, ceed_fallback;
109   CeedOperator op_fallback;
110 
111   // Check not already created
112   if (op->op_fallback) return CEED_ERROR_SUCCESS;
113 
114   // Fallback Ceed
115   CeedCall(CeedOperatorGetCeed(op, &ceed));
116   CeedCall(CeedGetOperatorFallbackCeed(ceed, &ceed_fallback));
117   CeedCall(CeedDestroy(&ceed));
118   if (!ceed_fallback) return CEED_ERROR_SUCCESS;
119 
120   CeedDebug256(CeedOperatorReturnCeed(op), 1, "---------- CeedOperator Fallback ----------\n");
121   CeedDebug(CeedOperatorReturnCeed(op), "Creating fallback CeedOperator\n");
122 
123   // Clone Op
124   CeedCall(CeedOperatorIsComposite(op, &is_composite));
125   if (is_composite) {
126     CeedInt       num_suboperators;
127     CeedOperator *sub_operators;
128 
129     CeedCall(CeedCompositeOperatorCreate(ceed_fallback, &op_fallback));
130     CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators));
131     CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
132     for (CeedInt i = 0; i < num_suboperators; i++) {
133       CeedOperator op_sub_fallback;
134 
135       CeedCall(CeedOperatorGetFallback(sub_operators[i], &op_sub_fallback));
136       CeedCall(CeedCompositeOperatorAddSub(op_fallback, op_sub_fallback));
137     }
138   } else {
139     bool               is_at_points = false;
140     CeedInt            num_input_fields, num_output_fields;
141     CeedQFunction      qf_fallback = NULL, dqf_fallback = NULL, dqfT_fallback = NULL;
142     CeedOperatorField *input_fields, *output_fields;
143 
144     CeedCall(CeedQFunctionCreateFallback(ceed_fallback, op->qf, &qf_fallback));
145     CeedCall(CeedQFunctionCreateFallback(ceed_fallback, op->dqf, &dqf_fallback));
146     CeedCall(CeedQFunctionCreateFallback(ceed_fallback, op->dqfT, &dqfT_fallback));
147     CeedCall(CeedOperatorIsAtPoints(op, &is_at_points));
148     if (is_at_points) {
149       CeedVector          points;
150       CeedElemRestriction rstr_points;
151 
152       CeedCall(CeedOperatorCreateAtPoints(ceed_fallback, qf_fallback, dqf_fallback, dqfT_fallback, &op_fallback));
153       CeedCall(CeedOperatorAtPointsGetPoints(op, &rstr_points, &points));
154       CeedCall(CeedOperatorAtPointsSetPoints(op_fallback, rstr_points, points));
155       CeedCall(CeedVectorDestroy(&points));
156       CeedCall(CeedElemRestrictionDestroy(&rstr_points));
157     } else {
158       CeedCall(CeedOperatorCreate(ceed_fallback, qf_fallback, dqf_fallback, dqfT_fallback, &op_fallback));
159     }
160     CeedCall(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
161     for (CeedInt i = 0; i < num_input_fields; i++) {
162       const char         *field_name;
163       CeedVector          vec;
164       CeedElemRestriction rstr;
165       CeedBasis           basis;
166 
167       CeedCall(CeedOperatorFieldGetData(input_fields[i], &field_name, &rstr, &basis, &vec));
168       CeedCall(CeedOperatorSetField(op_fallback, field_name, rstr, basis, vec));
169       CeedCall(CeedVectorDestroy(&vec));
170       CeedCall(CeedElemRestrictionDestroy(&rstr));
171       CeedCall(CeedBasisDestroy(&basis));
172     }
173     for (CeedInt i = 0; i < num_output_fields; i++) {
174       const char         *field_name;
175       CeedVector          vec;
176       CeedElemRestriction rstr;
177       CeedBasis           basis;
178 
179       CeedCall(CeedOperatorFieldGetData(output_fields[i], &field_name, &rstr, &basis, &vec));
180       CeedCall(CeedOperatorSetField(op_fallback, field_name, rstr, basis, vec));
181       CeedCall(CeedVectorDestroy(&vec));
182       CeedCall(CeedElemRestrictionDestroy(&rstr));
183       CeedCall(CeedBasisDestroy(&basis));
184     }
185     {
186       CeedQFunctionAssemblyData data;
187 
188       CeedCall(CeedOperatorGetQFunctionAssemblyData(op, &data));
189       CeedCall(CeedQFunctionAssemblyDataReferenceCopy(data, &op_fallback->qf_assembled));
190     }
191     // Cleanup
192     CeedCall(CeedQFunctionDestroy(&qf_fallback));
193     CeedCall(CeedQFunctionDestroy(&dqf_fallback));
194     CeedCall(CeedQFunctionDestroy(&dqfT_fallback));
195   }
196   CeedCall(CeedOperatorSetName(op_fallback, op->name));
197   CeedCall(CeedOperatorCheckReady(op_fallback));
198   // Note: No ref-counting here so we don't get caught in a reference loop.
199   //       The op holds the only reference to op_fallback and is responsible for deleting itself and op_fallback.
200   op->op_fallback                 = op_fallback;
201   op_fallback->op_fallback_parent = op;
202   CeedCall(CeedDestroy(&ceed_fallback));
203   return CEED_ERROR_SUCCESS;
204 }
205 
206 /**
207   @brief Core logic for assembling operator diagonal or point block diagonal
208 
209   @param[in]  op             `CeedOperator` to assemble diagonal or point block diagonal
210   @param[in]  request        Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
211   @param[in]  is_point_block Boolean flag to assemble diagonal or point block diagonal
212   @param[out] assembled      `CeedVector` to store assembled diagonal
213 
214   @return An error code: 0 - success, otherwise - failure
215 
216   @ref Developer
217 **/
218 static inline int CeedSingleOperatorLinearAssembleAddDiagonal_Mesh(CeedOperator op, CeedRequest *request, const bool is_point_block,
219                                                                    CeedVector assembled) {
220   bool is_composite;
221 
222   CeedCall(CeedOperatorIsComposite(op, &is_composite));
223   CeedCheck(!is_composite, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED, "Composite operator not supported");
224 
225   // Assemble QFunction
226   CeedInt             layout_qf[3];
227   const CeedScalar   *assembled_qf_array;
228   CeedVector          assembled_qf        = NULL;
229   CeedElemRestriction assembled_elem_rstr = NULL;
230 
231   CeedCall(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_elem_rstr, request));
232   CeedCall(CeedElemRestrictionGetELayout(assembled_elem_rstr, layout_qf));
233   CeedCall(CeedElemRestrictionDestroy(&assembled_elem_rstr));
234   CeedCall(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_HOST, &assembled_qf_array));
235 
236   // Get assembly data
237   const CeedEvalMode     **eval_modes_in, **eval_modes_out;
238   CeedInt                  num_active_bases_in, *num_eval_modes_in, num_active_bases_out, *num_eval_modes_out;
239   CeedSize               **eval_mode_offsets_in, **eval_mode_offsets_out, num_output_components;
240   CeedBasis               *active_bases_in, *active_bases_out;
241   CeedElemRestriction     *active_elem_rstrs_in, *active_elem_rstrs_out;
242   CeedOperatorAssemblyData data;
243 
244   CeedCall(CeedOperatorGetOperatorAssemblyData(op, &data));
245   CeedCall(CeedOperatorAssemblyDataGetEvalModes(data, &num_active_bases_in, &num_eval_modes_in, &eval_modes_in, &eval_mode_offsets_in,
246                                                 &num_active_bases_out, &num_eval_modes_out, &eval_modes_out, &eval_mode_offsets_out,
247                                                 &num_output_components));
248   CeedCall(CeedOperatorAssemblyDataGetBases(data, NULL, &active_bases_in, NULL, NULL, &active_bases_out, NULL));
249   CeedCall(CeedOperatorAssemblyDataGetElemRestrictions(data, NULL, &active_elem_rstrs_in, NULL, &active_elem_rstrs_out));
250 
251   // Loop over all active bases (find matching input/output pairs)
252   for (CeedInt b = 0; b < CeedIntMin(num_active_bases_in, num_active_bases_out); b++) {
253     CeedInt             b_in, b_out, num_elem, num_nodes, num_qpts, num_comp;
254     bool                has_eval_none = false;
255     CeedScalar         *elem_diag_array, *identity = NULL;
256     CeedVector          elem_diag;
257     CeedElemRestriction diag_elem_rstr;
258 
259     if (num_active_bases_in <= num_active_bases_out) {
260       b_in = b;
261       for (b_out = 0; b_out < num_active_bases_out; b_out++) {
262         if (active_bases_in[b_in] == active_bases_out[b_out]) {
263           break;
264         }
265       }
266       if (b_out == num_active_bases_out) {
267         continue;
268       }  // No matching output basis found
269     } else {
270       b_out = b;
271       for (b_in = 0; b_in < num_active_bases_in; b_in++) {
272         if (active_bases_in[b_in] == active_bases_out[b_out]) {
273           break;
274         }
275       }
276       if (b_in == num_active_bases_in) {
277         continue;
278       }  // No matching output basis found
279     }
280     CeedCheck(active_elem_rstrs_in[b_in] == active_elem_rstrs_out[b_out], CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED,
281               "Cannot assemble operator diagonal with different input and output active element restrictions");
282 
283     // Assemble point block diagonal restriction, if needed
284     if (is_point_block) {
285       CeedCall(CeedOperatorCreateActivePointBlockRestriction(active_elem_rstrs_in[b_in], &diag_elem_rstr));
286     } else {
287       CeedCall(CeedElemRestrictionCreateUnsignedCopy(active_elem_rstrs_in[b_in], &diag_elem_rstr));
288     }
289 
290     // Create diagonal vector
291     CeedCall(CeedElemRestrictionCreateVector(diag_elem_rstr, NULL, &elem_diag));
292 
293     // Assemble element operator diagonals
294     CeedCall(CeedVectorSetValue(elem_diag, 0.0));
295     CeedCall(CeedVectorGetArray(elem_diag, CEED_MEM_HOST, &elem_diag_array));
296     CeedCall(CeedElemRestrictionGetNumElements(diag_elem_rstr, &num_elem));
297     CeedCall(CeedBasisGetNumNodes(active_bases_in[b_in], &num_nodes));
298     CeedCall(CeedBasisGetNumComponents(active_bases_in[b_in], &num_comp));
299     if (active_bases_in[b_in] == CEED_BASIS_NONE) num_qpts = num_nodes;
300     else CeedCall(CeedBasisGetNumQuadraturePoints(active_bases_in[b_in], &num_qpts));
301 
302     // Construct identity matrix for basis if required
303     for (CeedInt i = 0; i < num_eval_modes_in[b_in]; i++) {
304       has_eval_none = has_eval_none || (eval_modes_in[b_in][i] == CEED_EVAL_NONE);
305     }
306     for (CeedInt i = 0; i < num_eval_modes_out[b_out]; i++) {
307       has_eval_none = has_eval_none || (eval_modes_out[b_out][i] == CEED_EVAL_NONE);
308     }
309     if (has_eval_none) {
310       CeedCall(CeedCalloc(num_qpts * num_nodes, &identity));
311       for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
312     }
313 
314     // Compute the diagonal of B^T D B
315     // Each element
316     for (CeedSize e = 0; e < num_elem; e++) {
317       // Each basis eval mode pair
318       CeedInt      d_out              = 0, q_comp_out;
319       CeedEvalMode eval_mode_out_prev = CEED_EVAL_NONE;
320 
321       for (CeedInt e_out = 0; e_out < num_eval_modes_out[b_out]; e_out++) {
322         CeedInt           d_in              = 0, q_comp_in;
323         const CeedScalar *B_t               = NULL;
324         CeedEvalMode      eval_mode_in_prev = CEED_EVAL_NONE;
325 
326         CeedCall(CeedOperatorGetBasisPointer(active_bases_out[b_out], eval_modes_out[b_out][e_out], identity, &B_t));
327         CeedCall(CeedBasisGetNumQuadratureComponents(active_bases_out[b_out], eval_modes_out[b_out][e_out], &q_comp_out));
328         if (q_comp_out > 1) {
329           if (e_out == 0 || eval_modes_out[b_out][e_out] != eval_mode_out_prev) d_out = 0;
330           else B_t = &B_t[(++d_out) * num_qpts * num_nodes];
331         }
332         eval_mode_out_prev = eval_modes_out[b_out][e_out];
333 
334         for (CeedInt e_in = 0; e_in < num_eval_modes_in[b_in]; e_in++) {
335           const CeedScalar *B = NULL;
336 
337           CeedCall(CeedOperatorGetBasisPointer(active_bases_in[b_in], eval_modes_in[b_in][e_in], identity, &B));
338           CeedCall(CeedBasisGetNumQuadratureComponents(active_bases_in[b_in], eval_modes_in[b_in][e_in], &q_comp_in));
339           if (q_comp_in > 1) {
340             if (e_in == 0 || eval_modes_in[b_in][e_in] != eval_mode_in_prev) d_in = 0;
341             else B = &B[(++d_in) * num_qpts * num_nodes];
342           }
343           eval_mode_in_prev = eval_modes_in[b_in][e_in];
344 
345           // Each component
346           for (CeedInt c_out = 0; c_out < num_comp; c_out++) {
347             // Each qpt/node pair
348             for (CeedInt q = 0; q < num_qpts; q++) {
349               if (is_point_block) {
350                 // Point Block Diagonal
351                 for (CeedInt c_in = 0; c_in < num_comp; c_in++) {
352                   const CeedSize c_offset =
353                       (eval_mode_offsets_in[b_in][e_in] + c_in) * num_output_components + eval_mode_offsets_out[b_out][e_out] + c_out;
354                   const CeedScalar qf_value = assembled_qf_array[q * layout_qf[0] + c_offset * layout_qf[1] + e * layout_qf[2]];
355 
356                   for (CeedInt n = 0; n < num_nodes; n++) {
357                     elem_diag_array[((e * num_comp + c_out) * num_comp + c_in) * num_nodes + n] +=
358                         B_t[q * num_nodes + n] * qf_value * B[q * num_nodes + n];
359                   }
360                 }
361               } else {
362                 // Diagonal Only
363                 const CeedInt c_offset =
364                     (eval_mode_offsets_in[b_in][e_in] + c_out) * num_output_components + eval_mode_offsets_out[b_out][e_out] + c_out;
365                 const CeedScalar qf_value = assembled_qf_array[q * layout_qf[0] + c_offset * layout_qf[1] + e * layout_qf[2]];
366 
367                 for (CeedInt n = 0; n < num_nodes; n++) {
368                   elem_diag_array[(e * num_comp + c_out) * num_nodes + n] += B_t[q * num_nodes + n] * qf_value * B[q * num_nodes + n];
369                 }
370               }
371             }
372           }
373         }
374       }
375     }
376     CeedCall(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
377 
378     // Assemble local operator diagonal
379     CeedCall(CeedElemRestrictionApply(diag_elem_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
380 
381     // Cleanup
382     CeedCall(CeedElemRestrictionDestroy(&diag_elem_rstr));
383     CeedCall(CeedVectorDestroy(&elem_diag));
384     CeedCall(CeedFree(&identity));
385   }
386   CeedCall(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
387   CeedCall(CeedVectorDestroy(&assembled_qf));
388   return CEED_ERROR_SUCCESS;
389 }
390 
391 /**
392   @brief Core logic for assembling operator diagonal or point block diagonal
393 
394   @param[in]  op             `CeedOperator` to assemble diagonal or point block diagonal
395   @param[in]  request        Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
396   @param[in]  is_point_block Boolean flag to assemble diagonal or point block diagonal
397   @param[out] assembled      `CeedVector` to store assembled diagonal
398 
399   @return An error code: 0 - success, otherwise - failure
400 
401   @ref Developer
402 **/
403 static inline int CeedSingleOperatorLinearAssembleAddDiagonal(CeedOperator op, CeedRequest *request, const bool is_point_block,
404                                                               CeedVector assembled) {
405   bool is_at_points;
406 
407   CeedCall(CeedOperatorIsAtPoints(op, &is_at_points));
408   CeedCheck(!is_at_points, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED, "AtPoints operator not supported");
409   CeedCall(CeedSingleOperatorLinearAssembleAddDiagonal_Mesh(op, request, is_point_block, assembled));
410   return CEED_ERROR_SUCCESS;
411 }
412 
413 /**
414   @brief Core logic for assembling composite operator diagonal
415 
416   @param[in]  op             `CeedOperator` to assemble point block diagonal
417   @param[in]  request        Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
418   @param[in]  is_point_block Boolean flag to assemble diagonal or point block diagonal
419   @param[out] assembled      `CeedVector` to store assembled diagonal
420 
421   @return An error code: 0 - success, otherwise - failure
422 
423   @ref Developer
424 **/
425 static inline int CeedCompositeOperatorLinearAssembleAddDiagonal(CeedOperator op, CeedRequest *request, const bool is_point_block,
426                                                                  CeedVector assembled) {
427   CeedInt       num_sub;
428   CeedOperator *suboperators;
429 
430   CeedCall(CeedCompositeOperatorGetNumSub(op, &num_sub));
431   CeedCall(CeedCompositeOperatorGetSubList(op, &suboperators));
432   for (CeedInt i = 0; i < num_sub; i++) {
433     if (is_point_block) {
434       CeedCall(CeedOperatorLinearAssembleAddPointBlockDiagonal(suboperators[i], assembled, request));
435     } else {
436       CeedCall(CeedOperatorLinearAssembleAddDiagonal(suboperators[i], assembled, request));
437     }
438   }
439   return CEED_ERROR_SUCCESS;
440 }
441 
442 /**
443   @brief Build nonzero pattern for non-composite CeedOperator`.
444 
445   Users should generally use @ref CeedOperatorLinearAssembleSymbolic().
446 
447   @param[in]  op     `CeedOperator` to assemble nonzero pattern
448   @param[in]  offset Offset for number of entries
449   @param[out] rows   Row number for each entry
450   @param[out] cols   Column number for each entry
451 
452   @return An error code: 0 - success, otherwise - failure
453 
454   @ref Developer
455 **/
456 static int CeedSingleOperatorAssembleSymbolic(CeedOperator op, CeedInt offset, CeedInt *rows, CeedInt *cols) {
457   Ceed                ceed;
458   bool                is_composite;
459   CeedSize            num_nodes_in, num_nodes_out, local_num_entries, count = 0;
460   CeedInt             num_elem_in, elem_size_in, num_comp_in, layout_er_in[3];
461   CeedInt             num_elem_out, elem_size_out, num_comp_out, layout_er_out[3];
462   CeedScalar         *array;
463   const CeedScalar   *elem_dof_a_in, *elem_dof_a_out;
464   CeedVector          index_vec_in, index_vec_out, elem_dof_in, elem_dof_out;
465   CeedElemRestriction elem_rstr_in, elem_rstr_out, index_elem_rstr_in, index_elem_rstr_out;
466 
467   CeedCall(CeedOperatorIsComposite(op, &is_composite));
468   CeedCall(CeedOperatorGetCeed(op, &ceed));
469   CeedCheck(!is_composite, ceed, CEED_ERROR_UNSUPPORTED, "Composite operator not supported");
470 
471   CeedCall(CeedOperatorGetActiveVectorLengths(op, &num_nodes_in, &num_nodes_out));
472   CeedCall(CeedOperatorGetActiveElemRestrictions(op, &elem_rstr_in, &elem_rstr_out));
473   CeedCall(CeedElemRestrictionGetNumElements(elem_rstr_in, &num_elem_in));
474   CeedCall(CeedElemRestrictionGetElementSize(elem_rstr_in, &elem_size_in));
475   CeedCall(CeedElemRestrictionGetNumComponents(elem_rstr_in, &num_comp_in));
476   CeedCall(CeedElemRestrictionGetELayout(elem_rstr_in, layout_er_in));
477 
478   // Determine elem_dof relation for input
479   CeedCall(CeedVectorCreate(ceed, num_nodes_in, &index_vec_in));
480   CeedCall(CeedVectorGetArrayWrite(index_vec_in, CEED_MEM_HOST, &array));
481   for (CeedSize i = 0; i < num_nodes_in; i++) array[i] = i;
482   CeedCall(CeedVectorRestoreArray(index_vec_in, &array));
483   CeedCall(CeedVectorCreate(ceed, num_elem_in * elem_size_in * num_comp_in, &elem_dof_in));
484   CeedCall(CeedVectorSetValue(elem_dof_in, 0.0));
485   CeedCall(CeedElemRestrictionCreateUnorientedCopy(elem_rstr_in, &index_elem_rstr_in));
486   CeedCall(CeedElemRestrictionApply(index_elem_rstr_in, CEED_NOTRANSPOSE, index_vec_in, elem_dof_in, CEED_REQUEST_IMMEDIATE));
487   CeedCall(CeedVectorGetArrayRead(elem_dof_in, CEED_MEM_HOST, &elem_dof_a_in));
488   CeedCall(CeedVectorDestroy(&index_vec_in));
489   CeedCall(CeedElemRestrictionDestroy(&index_elem_rstr_in));
490 
491   if (elem_rstr_in != elem_rstr_out) {
492     CeedCall(CeedElemRestrictionGetNumElements(elem_rstr_out, &num_elem_out));
493     CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED,
494               "Active input and output operator restrictions must have the same number of elements."
495               " Input has %" CeedInt_FMT " elements; output has %" CeedInt_FMT "elements.",
496               num_elem_in, num_elem_out);
497     CeedCall(CeedElemRestrictionGetElementSize(elem_rstr_out, &elem_size_out));
498     CeedCall(CeedElemRestrictionGetNumComponents(elem_rstr_out, &num_comp_out));
499     CeedCall(CeedElemRestrictionGetELayout(elem_rstr_out, layout_er_out));
500 
501     // Determine elem_dof relation for output
502     CeedCall(CeedVectorCreate(ceed, num_nodes_out, &index_vec_out));
503     CeedCall(CeedVectorGetArrayWrite(index_vec_out, CEED_MEM_HOST, &array));
504     for (CeedSize i = 0; i < num_nodes_out; i++) array[i] = i;
505     CeedCall(CeedVectorRestoreArray(index_vec_out, &array));
506     CeedCall(CeedVectorCreate(ceed, num_elem_out * elem_size_out * num_comp_out, &elem_dof_out));
507     CeedCall(CeedVectorSetValue(elem_dof_out, 0.0));
508     CeedCall(CeedElemRestrictionCreateUnorientedCopy(elem_rstr_out, &index_elem_rstr_out));
509     CeedCall(CeedElemRestrictionApply(index_elem_rstr_out, CEED_NOTRANSPOSE, index_vec_out, elem_dof_out, CEED_REQUEST_IMMEDIATE));
510     CeedCall(CeedVectorGetArrayRead(elem_dof_out, CEED_MEM_HOST, &elem_dof_a_out));
511     CeedCall(CeedVectorDestroy(&index_vec_out));
512     CeedCall(CeedElemRestrictionDestroy(&index_elem_rstr_out));
513   } else {
514     num_elem_out     = num_elem_in;
515     elem_size_out    = elem_size_in;
516     num_comp_out     = num_comp_in;
517     layout_er_out[0] = layout_er_in[0];
518     layout_er_out[1] = layout_er_in[1];
519     layout_er_out[2] = layout_er_in[2];
520     elem_dof_a_out   = elem_dof_a_in;
521   }
522   local_num_entries = (CeedSize)elem_size_out * num_comp_out * elem_size_in * num_comp_in * num_elem_in;
523 
524   // Determine i, j locations for element matrices
525   for (CeedInt e = 0; e < num_elem_in; e++) {
526     for (CeedInt comp_in = 0; comp_in < num_comp_in; comp_in++) {
527       for (CeedInt comp_out = 0; comp_out < num_comp_out; comp_out++) {
528         for (CeedInt i = 0; i < elem_size_out; i++) {
529           for (CeedInt j = 0; j < elem_size_in; j++) {
530             const CeedInt elem_dof_index_row = i * layout_er_out[0] + comp_out * layout_er_out[1] + e * layout_er_out[2];
531             const CeedInt elem_dof_index_col = j * layout_er_in[0] + comp_in * layout_er_in[1] + e * layout_er_in[2];
532             const CeedInt row                = elem_dof_a_out[elem_dof_index_row];
533             const CeedInt col                = elem_dof_a_in[elem_dof_index_col];
534 
535             rows[offset + count] = row;
536             cols[offset + count] = col;
537             count++;
538           }
539         }
540       }
541     }
542   }
543   CeedCheck(count == local_num_entries, ceed, CEED_ERROR_MAJOR, "Error computing assembled entries");
544   CeedCall(CeedVectorRestoreArrayRead(elem_dof_in, &elem_dof_a_in));
545   CeedCall(CeedVectorDestroy(&elem_dof_in));
546   if (elem_rstr_in != elem_rstr_out) {
547     CeedCall(CeedVectorRestoreArrayRead(elem_dof_out, &elem_dof_a_out));
548     CeedCall(CeedVectorDestroy(&elem_dof_out));
549   }
550   CeedCall(CeedElemRestrictionDestroy(&elem_rstr_in));
551   CeedCall(CeedElemRestrictionDestroy(&elem_rstr_out));
552   CeedCall(CeedDestroy(&ceed));
553   return CEED_ERROR_SUCCESS;
554 }
555 
556 /**
557   @brief Assemble nonzero entries for non-composite `CeedOperator`.
558 
559   Users should generally use @ref CeedOperatorLinearAssemble().
560 
561   @param[in]  op     `CeedOperator` to assemble
562   @param[in]  offset Offset for number of entries
563   @param[out] values Values to assemble into matrix
564 
565   @return An error code: 0 - success, otherwise - failure
566 
567   @ref Developer
568 **/
569 static int CeedSingleOperatorAssemble(CeedOperator op, CeedInt offset, CeedVector values) {
570   bool is_composite;
571 
572   CeedCall(CeedOperatorIsComposite(op, &is_composite));
573   CeedCheck(!is_composite, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED, "Composite operator not supported");
574 
575   // Early exit for empty operator
576   {
577     CeedInt num_elem = 0;
578 
579     CeedCall(CeedOperatorGetNumElements(op, &num_elem));
580     if (num_elem == 0) return CEED_ERROR_SUCCESS;
581   }
582 
583   if (op->LinearAssembleSingle) {
584     // Backend version
585     CeedCall(op->LinearAssembleSingle(op, offset, values));
586     return CEED_ERROR_SUCCESS;
587   } else {
588     // Operator fallback
589     CeedOperator op_fallback;
590 
591     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
592     if (op_fallback) {
593       CeedCall(CeedSingleOperatorAssemble(op_fallback, offset, values));
594       return CEED_ERROR_SUCCESS;
595     }
596   }
597 
598   // Assemble QFunction
599   CeedInt             layout_qf[3];
600   const CeedScalar   *assembled_qf_array;
601   CeedVector          assembled_qf        = NULL;
602   CeedElemRestriction assembled_elem_rstr = NULL;
603 
604   CeedCall(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_elem_rstr, CEED_REQUEST_IMMEDIATE));
605   CeedCall(CeedElemRestrictionGetELayout(assembled_elem_rstr, layout_qf));
606   CeedCall(CeedElemRestrictionDestroy(&assembled_elem_rstr));
607   CeedCall(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_HOST, &assembled_qf_array));
608 
609   // Get assembly data
610   CeedInt                  num_elem_in, elem_size_in, num_comp_in, num_qpts_in;
611   CeedInt                  num_elem_out, elem_size_out, num_comp_out, num_qpts_out;
612   CeedSize                 local_num_entries, count = 0;
613   const CeedEvalMode     **eval_modes_in, **eval_modes_out;
614   CeedInt                  num_active_bases_in, *num_eval_modes_in, num_active_bases_out, *num_eval_modes_out;
615   CeedBasis               *active_bases_in, *active_bases_out, basis_in, basis_out;
616   const CeedScalar       **B_mats_in, **B_mats_out, *B_mat_in, *B_mat_out;
617   CeedElemRestriction      elem_rstr_in, elem_rstr_out;
618   CeedRestrictionType      elem_rstr_type_in, elem_rstr_type_out;
619   const bool              *elem_rstr_orients_in = NULL, *elem_rstr_orients_out = NULL;
620   const CeedInt8          *elem_rstr_curl_orients_in = NULL, *elem_rstr_curl_orients_out = NULL;
621   CeedOperatorAssemblyData data;
622 
623   CeedCall(CeedOperatorGetOperatorAssemblyData(op, &data));
624   CeedCall(CeedOperatorAssemblyDataGetEvalModes(data, &num_active_bases_in, &num_eval_modes_in, &eval_modes_in, NULL, &num_active_bases_out,
625                                                 &num_eval_modes_out, &eval_modes_out, NULL, NULL));
626 
627   CeedCheck(num_active_bases_in == 1 && num_active_bases_out == 1, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED,
628             "Cannot assemble operator with multiple active bases");
629   CeedCheck(num_eval_modes_in[0] > 0 && num_eval_modes_out[0] > 0, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED,
630             "Cannot assemble operator without inputs/outputs");
631 
632   CeedCall(CeedOperatorAssemblyDataGetBases(data, NULL, &active_bases_in, &B_mats_in, NULL, &active_bases_out, &B_mats_out));
633   CeedCall(CeedOperatorGetActiveElemRestrictions(op, &elem_rstr_in, &elem_rstr_out));
634   basis_in  = active_bases_in[0];
635   basis_out = active_bases_out[0];
636   B_mat_in  = B_mats_in[0];
637   B_mat_out = B_mats_out[0];
638 
639   CeedCall(CeedElemRestrictionGetNumElements(elem_rstr_in, &num_elem_in));
640   CeedCall(CeedElemRestrictionGetElementSize(elem_rstr_in, &elem_size_in));
641   CeedCall(CeedElemRestrictionGetNumComponents(elem_rstr_in, &num_comp_in));
642   if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in;
643   else CeedCall(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in));
644 
645   CeedCall(CeedElemRestrictionGetType(elem_rstr_in, &elem_rstr_type_in));
646   if (elem_rstr_type_in == CEED_RESTRICTION_ORIENTED) {
647     CeedCall(CeedElemRestrictionGetOrientations(elem_rstr_in, CEED_MEM_HOST, &elem_rstr_orients_in));
648   } else if (elem_rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
649     CeedCall(CeedElemRestrictionGetCurlOrientations(elem_rstr_in, CEED_MEM_HOST, &elem_rstr_curl_orients_in));
650   }
651 
652   if (elem_rstr_in != elem_rstr_out) {
653     CeedCall(CeedElemRestrictionGetNumElements(elem_rstr_out, &num_elem_out));
654     CeedCheck(num_elem_in == num_elem_out, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED,
655               "Active input and output operator restrictions must have the same number of elements."
656               " Input has %" CeedInt_FMT " elements; output has %" CeedInt_FMT "elements.",
657               num_elem_in, num_elem_out);
658     CeedCall(CeedElemRestrictionGetElementSize(elem_rstr_out, &elem_size_out));
659     CeedCall(CeedElemRestrictionGetNumComponents(elem_rstr_out, &num_comp_out));
660     if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out;
661     else CeedCall(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out));
662     CeedCheck(num_qpts_in == num_qpts_out, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED,
663               "Active input and output bases must have the same number of quadrature points."
664               " Input has %" CeedInt_FMT " points; output has %" CeedInt_FMT "points.",
665               num_qpts_in, num_qpts_out);
666 
667     CeedCall(CeedElemRestrictionGetType(elem_rstr_out, &elem_rstr_type_out));
668     if (elem_rstr_type_out == CEED_RESTRICTION_ORIENTED) {
669       CeedCall(CeedElemRestrictionGetOrientations(elem_rstr_out, CEED_MEM_HOST, &elem_rstr_orients_out));
670     } else if (elem_rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
671       CeedCall(CeedElemRestrictionGetCurlOrientations(elem_rstr_out, CEED_MEM_HOST, &elem_rstr_curl_orients_out));
672     }
673   } else {
674     num_elem_out  = num_elem_in;
675     elem_size_out = elem_size_in;
676     num_comp_out  = num_comp_in;
677     num_qpts_out  = num_qpts_in;
678 
679     elem_rstr_orients_out      = elem_rstr_orients_in;
680     elem_rstr_curl_orients_out = elem_rstr_curl_orients_in;
681   }
682   local_num_entries = (CeedSize)elem_size_out * num_comp_out * elem_size_in * num_comp_in * num_elem_in;
683 
684   // Loop over elements and put in data structure
685   // We store B_mat_in, B_mat_out, BTD, elem_mat in row-major order
686   CeedTensorContract contract;
687   CeedScalar        *vals, *BTD_mat = NULL, *elem_mat = NULL, *elem_mat_b = NULL;
688 
689   CeedCall(CeedBasisGetTensorContract(basis_in, &contract));
690   CeedCall(CeedCalloc(elem_size_out * num_qpts_in * num_eval_modes_in[0], &BTD_mat));
691   CeedCall(CeedCalloc(elem_size_out * elem_size_in, &elem_mat));
692   if (elem_rstr_curl_orients_in || elem_rstr_curl_orients_out) CeedCall(CeedCalloc(elem_size_out * elem_size_in, &elem_mat_b));
693 
694   CeedCall(CeedVectorGetArray(values, CEED_MEM_HOST, &vals));
695   for (CeedSize e = 0; e < num_elem_in; e++) {
696     for (CeedInt comp_in = 0; comp_in < num_comp_in; comp_in++) {
697       for (CeedInt comp_out = 0; comp_out < num_comp_out; comp_out++) {
698         // Compute B^T*D
699         for (CeedSize n = 0; n < elem_size_out; n++) {
700           for (CeedSize q = 0; q < num_qpts_in; q++) {
701             for (CeedInt e_in = 0; e_in < num_eval_modes_in[0]; e_in++) {
702               const CeedSize btd_index = n * (num_qpts_in * num_eval_modes_in[0]) + q * num_eval_modes_in[0] + e_in;
703               CeedScalar     sum       = 0.0;
704 
705               for (CeedInt e_out = 0; e_out < num_eval_modes_out[0]; e_out++) {
706                 const CeedSize b_out_index     = (q * num_eval_modes_out[0] + e_out) * elem_size_out + n;
707                 const CeedSize eval_mode_index = ((e_in * num_comp_in + comp_in) * num_eval_modes_out[0] + e_out) * num_comp_out + comp_out;
708                 const CeedSize qf_index        = q * layout_qf[0] + eval_mode_index * layout_qf[1] + e * layout_qf[2];
709 
710                 sum += B_mat_out[b_out_index] * assembled_qf_array[qf_index];
711               }
712               BTD_mat[btd_index] = sum;
713             }
714           }
715         }
716 
717         // Form element matrix itself (for each block component)
718         if (contract) {
719           CeedCall(CeedTensorContractApply(contract, 1, num_qpts_in * num_eval_modes_in[0], elem_size_in, elem_size_out, BTD_mat, CEED_NOTRANSPOSE,
720                                            false, B_mat_in, elem_mat));
721         } else {
722           Ceed ceed;
723 
724           CeedCall(CeedOperatorGetCeed(op, &ceed));
725           CeedCall(CeedMatrixMatrixMultiply(ceed, BTD_mat, B_mat_in, elem_mat, elem_size_out, elem_size_in, num_qpts_in * num_eval_modes_in[0]));
726           CeedCall(CeedDestroy(&ceed));
727         }
728 
729         // Transform the element matrix if required
730         if (elem_rstr_orients_out) {
731           const bool *elem_orients = &elem_rstr_orients_out[e * elem_size_out];
732 
733           for (CeedInt i = 0; i < elem_size_out; i++) {
734             const double orient = elem_orients[i] ? -1.0 : 1.0;
735 
736             for (CeedInt j = 0; j < elem_size_in; j++) {
737               elem_mat[i * elem_size_in + j] *= orient;
738             }
739           }
740         } else if (elem_rstr_curl_orients_out) {
741           const CeedInt8 *elem_curl_orients = &elem_rstr_curl_orients_out[e * 3 * elem_size_out];
742 
743           // T^T*(B^T*D*B)
744           memcpy(elem_mat_b, elem_mat, elem_size_out * elem_size_in * sizeof(CeedScalar));
745           for (CeedInt i = 0; i < elem_size_out; i++) {
746             for (CeedInt j = 0; j < elem_size_in; j++) {
747               elem_mat[i * elem_size_in + j] = elem_mat_b[i * elem_size_in + j] * elem_curl_orients[3 * i + 1] +
748                                                (i > 0 ? elem_mat_b[(i - 1) * elem_size_in + j] * elem_curl_orients[3 * i - 1] : 0.0) +
749                                                (i < elem_size_out - 1 ? elem_mat_b[(i + 1) * elem_size_in + j] * elem_curl_orients[3 * i + 3] : 0.0);
750             }
751           }
752         }
753         if (elem_rstr_orients_in) {
754           const bool *elem_orients = &elem_rstr_orients_in[e * elem_size_in];
755 
756           for (CeedInt i = 0; i < elem_size_out; i++) {
757             for (CeedInt j = 0; j < elem_size_in; j++) {
758               elem_mat[i * elem_size_in + j] *= elem_orients[j] ? -1.0 : 1.0;
759             }
760           }
761         } else if (elem_rstr_curl_orients_in) {
762           const CeedInt8 *elem_curl_orients = &elem_rstr_curl_orients_in[e * 3 * elem_size_in];
763 
764           // (B^T*D*B)*T
765           memcpy(elem_mat_b, elem_mat, elem_size_out * elem_size_in * sizeof(CeedScalar));
766           for (CeedInt i = 0; i < elem_size_out; i++) {
767             for (CeedInt j = 0; j < elem_size_in; j++) {
768               elem_mat[i * elem_size_in + j] = elem_mat_b[i * elem_size_in + j] * elem_curl_orients[3 * j + 1] +
769                                                (j > 0 ? elem_mat_b[i * elem_size_in + j - 1] * elem_curl_orients[3 * j - 1] : 0.0) +
770                                                (j < elem_size_in - 1 ? elem_mat_b[i * elem_size_in + j + 1] * elem_curl_orients[3 * j + 3] : 0.0);
771             }
772           }
773         }
774 
775         // Put element matrix in coordinate data structure
776         for (CeedInt i = 0; i < elem_size_out; i++) {
777           for (CeedInt j = 0; j < elem_size_in; j++) {
778             vals[offset + count] = elem_mat[i * elem_size_in + j];
779             count++;
780           }
781         }
782       }
783     }
784   }
785   CeedCheck(count == local_num_entries, CeedOperatorReturnCeed(op), CEED_ERROR_MAJOR, "Error computing entries");
786   CeedCall(CeedVectorRestoreArray(values, &vals));
787 
788   // Cleanup
789   CeedCall(CeedFree(&BTD_mat));
790   CeedCall(CeedFree(&elem_mat));
791   CeedCall(CeedFree(&elem_mat_b));
792   if (elem_rstr_type_in == CEED_RESTRICTION_ORIENTED) {
793     CeedCall(CeedElemRestrictionRestoreOrientations(elem_rstr_in, &elem_rstr_orients_in));
794   } else if (elem_rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
795     CeedCall(CeedElemRestrictionRestoreCurlOrientations(elem_rstr_in, &elem_rstr_curl_orients_in));
796   }
797   if (elem_rstr_in != elem_rstr_out) {
798     if (elem_rstr_type_out == CEED_RESTRICTION_ORIENTED) {
799       CeedCall(CeedElemRestrictionRestoreOrientations(elem_rstr_out, &elem_rstr_orients_out));
800     } else if (elem_rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
801       CeedCall(CeedElemRestrictionRestoreCurlOrientations(elem_rstr_out, &elem_rstr_curl_orients_out));
802     }
803   }
804   CeedCall(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
805   CeedCall(CeedVectorDestroy(&assembled_qf));
806   CeedCall(CeedElemRestrictionDestroy(&elem_rstr_in));
807   CeedCall(CeedElemRestrictionDestroy(&elem_rstr_out));
808   return CEED_ERROR_SUCCESS;
809 }
810 
811 /**
812   @brief Count number of entries for assembled `CeedOperator`
813 
814   @param[in]  op          `CeedOperator` to assemble
815   @param[out] num_entries Number of entries in assembled representation
816 
817   @return An error code: 0 - success, otherwise - failure
818 
819   @ref Utility
820 **/
821 static int CeedSingleOperatorAssemblyCountEntries(CeedOperator op, CeedSize *num_entries) {
822   bool                is_composite;
823   CeedInt             num_elem_in, elem_size_in, num_comp_in, num_elem_out, elem_size_out, num_comp_out;
824   CeedElemRestriction rstr_in, rstr_out;
825 
826   CeedCall(CeedOperatorIsComposite(op, &is_composite));
827   CeedCheck(!is_composite, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED, "Composite operator not supported");
828 
829   CeedCall(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
830   CeedCall(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in));
831   CeedCall(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
832   CeedCall(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in));
833   if (rstr_in != rstr_out) {
834     CeedCall(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out));
835     CeedCheck(num_elem_in == num_elem_out, CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED,
836               "Active input and output operator restrictions must have the same number of elements."
837               " Input has %" CeedInt_FMT " elements; output has %" CeedInt_FMT "elements.",
838               num_elem_in, num_elem_out);
839     CeedCall(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
840     CeedCall(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out));
841   } else {
842     num_elem_out  = num_elem_in;
843     elem_size_out = elem_size_in;
844     num_comp_out  = num_comp_in;
845   }
846   CeedCall(CeedElemRestrictionDestroy(&rstr_in));
847   CeedCall(CeedElemRestrictionDestroy(&rstr_out));
848   *num_entries = (CeedSize)elem_size_in * num_comp_in * elem_size_out * num_comp_out * num_elem_in;
849   return CEED_ERROR_SUCCESS;
850 }
851 
852 /**
853   @brief Common code for creating a multigrid coarse `CeedOperator` and level transfer `CeedOperator` for a `CeedOperator`
854 
855   @param[in]  op_fine      Fine grid `CeedOperator`
856   @param[in]  p_mult_fine  L-vector multiplicity in parallel gather/scatter, or `NULL` if not creating prolongation/restriction `CeedOperator`
857   @param[in]  rstr_coarse  Coarse grid `CeedElemRestriction`
858   @param[in]  basis_coarse Coarse grid active vector `CeedBasis`
859   @param[in]  basis_c_to_f `CeedBasis` for coarse to fine interpolation, or `NULL` if not creating prolongation/restriction operators
860   @param[out] op_coarse    Coarse grid `CeedOperator`
861   @param[out] op_prolong   Coarse to fine `CeedOperator`, or `NULL`
862   @param[out] op_restrict  Fine to coarse `CeedOperator`, or `NULL`
863 
864   @return An error code: 0 - success, otherwise - failure
865 
866   @ref Developer
867 **/
868 static int CeedSingleOperatorMultigridLevel(CeedOperator op_fine, CeedVector p_mult_fine, CeedElemRestriction rstr_coarse, CeedBasis basis_coarse,
869                                             CeedBasis basis_c_to_f, CeedOperator *op_coarse, CeedOperator *op_prolong, CeedOperator *op_restrict) {
870   bool                is_composite;
871   Ceed                ceed;
872   CeedInt             num_comp, num_input_fields, num_output_fields;
873   CeedVector          mult_vec         = NULL;
874   CeedElemRestriction rstr_p_mult_fine = NULL, rstr_fine = NULL;
875   CeedOperatorField  *input_fields, *output_fields;
876 
877   CeedCall(CeedOperatorGetCeed(op_fine, &ceed));
878 
879   // Check for composite operator
880   CeedCall(CeedOperatorIsComposite(op_fine, &is_composite));
881   CeedCheck(!is_composite, ceed, CEED_ERROR_UNSUPPORTED, "Automatic multigrid setup for composite operators not supported");
882 
883   // Coarse Grid
884   CeedCall(CeedOperatorCreate(ceed, op_fine->qf, op_fine->dqf, op_fine->dqfT, op_coarse));
885   CeedCall(CeedOperatorGetFields(op_fine, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
886   // -- Clone input fields
887   for (CeedInt i = 0; i < num_input_fields; i++) {
888     const char         *field_name;
889     CeedVector          vec;
890     CeedElemRestriction rstr  = NULL;
891     CeedBasis           basis = NULL;
892 
893     CeedCall(CeedOperatorFieldGetName(input_fields[i], &field_name));
894     CeedCall(CeedOperatorFieldGetVector(input_fields[i], &vec));
895     if (vec == CEED_VECTOR_ACTIVE) {
896       CeedCall(CeedElemRestrictionReferenceCopy(rstr_coarse, &rstr));
897       CeedCall(CeedBasisReferenceCopy(basis_coarse, &basis));
898       if (!rstr_fine) CeedCall(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_fine));
899     } else {
900       CeedCall(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr));
901       CeedCall(CeedOperatorFieldGetBasis(input_fields[i], &basis));
902     }
903     CeedCall(CeedOperatorSetField(*op_coarse, field_name, rstr, basis, vec));
904     CeedCall(CeedVectorDestroy(&vec));
905     CeedCall(CeedElemRestrictionDestroy(&rstr));
906     CeedCall(CeedBasisDestroy(&basis));
907   }
908   // -- Clone output fields
909   for (CeedInt i = 0; i < num_output_fields; i++) {
910     const char         *field_name;
911     CeedVector          vec;
912     CeedElemRestriction rstr  = NULL;
913     CeedBasis           basis = NULL;
914 
915     CeedCall(CeedOperatorFieldGetName(output_fields[i], &field_name));
916     CeedCall(CeedOperatorFieldGetVector(output_fields[i], &vec));
917     if (vec == CEED_VECTOR_ACTIVE) {
918       CeedCall(CeedElemRestrictionReferenceCopy(rstr_coarse, &rstr));
919       CeedCall(CeedBasisReferenceCopy(basis_coarse, &basis));
920       if (!rstr_fine) CeedCall(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_fine));
921     } else {
922       CeedCall(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr));
923       CeedCall(CeedOperatorFieldGetBasis(output_fields[i], &basis));
924     }
925     CeedCall(CeedOperatorSetField(*op_coarse, field_name, rstr, basis, vec));
926     CeedCall(CeedVectorDestroy(&vec));
927     CeedCall(CeedElemRestrictionDestroy(&rstr));
928     CeedCall(CeedBasisDestroy(&basis));
929   }
930   // -- Clone QFunctionAssemblyData
931   {
932     CeedQFunctionAssemblyData fine_data;
933 
934     CeedCall(CeedOperatorGetQFunctionAssemblyData(op_fine, &fine_data));
935     CeedCall(CeedQFunctionAssemblyDataReferenceCopy(fine_data, &(*op_coarse)->qf_assembled));
936   }
937 
938   // Multiplicity vector
939   if (op_restrict || op_prolong) {
940     CeedVector          mult_e_vec;
941     CeedRestrictionType rstr_type;
942 
943     CeedCall(CeedElemRestrictionGetType(rstr_fine, &rstr_type));
944     CeedCheck(rstr_type != CEED_RESTRICTION_CURL_ORIENTED, ceed, CEED_ERROR_UNSUPPORTED,
945               "Element restrictions created with CeedElemRestrictionCreateCurlOriented are not supported");
946     CeedCheck(p_mult_fine, ceed, CEED_ERROR_INCOMPATIBLE, "Prolongation or restriction operator creation requires fine grid multiplicity vector");
947     CeedCall(CeedElemRestrictionCreateUnsignedCopy(rstr_fine, &rstr_p_mult_fine));
948     CeedCall(CeedElemRestrictionCreateVector(rstr_fine, &mult_vec, &mult_e_vec));
949     CeedCall(CeedVectorSetValue(mult_e_vec, 0.0));
950     CeedCall(CeedElemRestrictionApply(rstr_p_mult_fine, CEED_NOTRANSPOSE, p_mult_fine, mult_e_vec, CEED_REQUEST_IMMEDIATE));
951     CeedCall(CeedVectorSetValue(mult_vec, 0.0));
952     CeedCall(CeedElemRestrictionApply(rstr_p_mult_fine, CEED_TRANSPOSE, mult_e_vec, mult_vec, CEED_REQUEST_IMMEDIATE));
953     CeedCall(CeedVectorDestroy(&mult_e_vec));
954     CeedCall(CeedVectorReciprocal(mult_vec));
955   }
956 
957   // Clone name
958   bool   has_name = op_fine->name;
959   size_t name_len = op_fine->name ? strlen(op_fine->name) : 0;
960   CeedCall(CeedOperatorSetName(*op_coarse, op_fine->name));
961 
962   // Check that coarse to fine basis is provided if prolong/restrict operators are requested
963   CeedCheck(basis_c_to_f || (!op_restrict && !op_prolong), ceed, CEED_ERROR_INCOMPATIBLE,
964             "Prolongation or restriction operator creation requires coarse-to-fine basis");
965 
966   // Restriction/Prolongation Operators
967   CeedCall(CeedBasisGetNumComponents(basis_coarse, &num_comp));
968 
969   // Restriction
970   if (op_restrict) {
971     CeedInt             *num_comp_r_data;
972     CeedQFunctionContext ctx_r;
973     CeedQFunction        qf_restrict;
974 
975     CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Scale", &qf_restrict));
976     CeedCall(CeedCalloc(1, &num_comp_r_data));
977     num_comp_r_data[0] = num_comp;
978     CeedCall(CeedQFunctionContextCreate(ceed, &ctx_r));
979     CeedCall(CeedQFunctionContextSetData(ctx_r, CEED_MEM_HOST, CEED_OWN_POINTER, sizeof(*num_comp_r_data), num_comp_r_data));
980     CeedCall(CeedQFunctionSetContext(qf_restrict, ctx_r));
981     CeedCall(CeedQFunctionContextDestroy(&ctx_r));
982     CeedCall(CeedQFunctionAddInput(qf_restrict, "input", num_comp, CEED_EVAL_NONE));
983     CeedCall(CeedQFunctionAddInput(qf_restrict, "scale", num_comp, CEED_EVAL_NONE));
984     CeedCall(CeedQFunctionAddOutput(qf_restrict, "output", num_comp, CEED_EVAL_INTERP));
985     CeedCall(CeedQFunctionSetUserFlopsEstimate(qf_restrict, num_comp));
986 
987     CeedCall(CeedOperatorCreate(ceed, qf_restrict, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, op_restrict));
988     CeedCall(CeedOperatorSetField(*op_restrict, "input", rstr_fine, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
989     CeedCall(CeedOperatorSetField(*op_restrict, "scale", rstr_p_mult_fine, CEED_BASIS_NONE, mult_vec));
990     CeedCall(CeedOperatorSetField(*op_restrict, "output", rstr_coarse, basis_c_to_f, CEED_VECTOR_ACTIVE));
991 
992     // Set name
993     char *restriction_name;
994 
995     CeedCall(CeedCalloc(17 + name_len, &restriction_name));
996     sprintf(restriction_name, "restriction%s%s", has_name ? " for " : "", has_name ? op_fine->name : "");
997     CeedCall(CeedOperatorSetName(*op_restrict, restriction_name));
998     CeedCall(CeedFree(&restriction_name));
999 
1000     // Check
1001     CeedCall(CeedOperatorCheckReady(*op_restrict));
1002 
1003     // Cleanup
1004     CeedCall(CeedQFunctionDestroy(&qf_restrict));
1005   }
1006 
1007   // Prolongation
1008   if (op_prolong) {
1009     CeedInt             *num_comp_p_data;
1010     CeedQFunctionContext ctx_p;
1011     CeedQFunction        qf_prolong;
1012 
1013     CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Scale", &qf_prolong));
1014     CeedCall(CeedCalloc(1, &num_comp_p_data));
1015     num_comp_p_data[0] = num_comp;
1016     CeedCall(CeedQFunctionContextCreate(ceed, &ctx_p));
1017     CeedCall(CeedQFunctionContextSetData(ctx_p, CEED_MEM_HOST, CEED_OWN_POINTER, sizeof(*num_comp_p_data), num_comp_p_data));
1018     CeedCall(CeedQFunctionSetContext(qf_prolong, ctx_p));
1019     CeedCall(CeedQFunctionContextDestroy(&ctx_p));
1020     CeedCall(CeedQFunctionAddInput(qf_prolong, "input", num_comp, CEED_EVAL_INTERP));
1021     CeedCall(CeedQFunctionAddInput(qf_prolong, "scale", num_comp, CEED_EVAL_NONE));
1022     CeedCall(CeedQFunctionAddOutput(qf_prolong, "output", num_comp, CEED_EVAL_NONE));
1023     CeedCall(CeedQFunctionSetUserFlopsEstimate(qf_prolong, num_comp));
1024 
1025     CeedCall(CeedOperatorCreate(ceed, qf_prolong, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, op_prolong));
1026     CeedCall(CeedOperatorSetField(*op_prolong, "input", rstr_coarse, basis_c_to_f, CEED_VECTOR_ACTIVE));
1027     CeedCall(CeedOperatorSetField(*op_prolong, "scale", rstr_p_mult_fine, CEED_BASIS_NONE, mult_vec));
1028     CeedCall(CeedOperatorSetField(*op_prolong, "output", rstr_fine, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1029 
1030     // Set name
1031     char *prolongation_name;
1032 
1033     CeedCall(CeedCalloc(18 + name_len, &prolongation_name));
1034     sprintf(prolongation_name, "prolongation%s%s", has_name ? " for " : "", has_name ? op_fine->name : "");
1035     CeedCall(CeedOperatorSetName(*op_prolong, prolongation_name));
1036     CeedCall(CeedFree(&prolongation_name));
1037 
1038     // Check
1039     CeedCall(CeedOperatorCheckReady(*op_prolong));
1040 
1041     // Cleanup
1042     CeedCall(CeedQFunctionDestroy(&qf_prolong));
1043   }
1044 
1045   // Check
1046   CeedCall(CeedOperatorCheckReady(*op_coarse));
1047 
1048   // Cleanup
1049   CeedCall(CeedDestroy(&ceed));
1050   CeedCall(CeedVectorDestroy(&mult_vec));
1051   CeedCall(CeedElemRestrictionDestroy(&rstr_fine));
1052   CeedCall(CeedElemRestrictionDestroy(&rstr_p_mult_fine));
1053   CeedCall(CeedBasisDestroy(&basis_c_to_f));
1054   return CEED_ERROR_SUCCESS;
1055 }
1056 
1057 /**
1058   @brief Build 1D mass matrix and Laplacian with perturbation
1059 
1060   @param[in]  interp_1d   Interpolation matrix in one dimension
1061   @param[in]  grad_1d     Gradient matrix in one dimension
1062   @param[in]  q_weight_1d Quadrature weights in one dimension
1063   @param[in]  P_1d        Number of basis nodes in one dimension
1064   @param[in]  Q_1d        Number of quadrature points in one dimension
1065   @param[in]  dim         Dimension of basis
1066   @param[out] mass        Assembled mass matrix in one dimension
1067   @param[out] laplace     Assembled perturbed Laplacian in one dimension
1068 
1069   @return An error code: 0 - success, otherwise - failure
1070 
1071   @ref Developer
1072 **/
1073 CeedPragmaOptimizeOff
1074 static int CeedBuildMassLaplace(const CeedScalar *interp_1d, const CeedScalar *grad_1d, const CeedScalar *q_weight_1d, CeedInt P_1d, CeedInt Q_1d,
1075                                 CeedInt dim, CeedScalar *mass, CeedScalar *laplace) {
1076   for (CeedInt i = 0; i < P_1d; i++) {
1077     for (CeedInt j = 0; j < P_1d; j++) {
1078       CeedScalar sum = 0.0;
1079       for (CeedInt k = 0; k < Q_1d; k++) sum += interp_1d[k * P_1d + i] * q_weight_1d[k] * interp_1d[k * P_1d + j];
1080       mass[i + j * P_1d] = sum;
1081     }
1082   }
1083   // -- Laplacian
1084   for (CeedInt i = 0; i < P_1d; i++) {
1085     for (CeedInt j = 0; j < P_1d; j++) {
1086       CeedScalar sum = 0.0;
1087 
1088       for (CeedInt k = 0; k < Q_1d; k++) sum += grad_1d[k * P_1d + i] * q_weight_1d[k] * grad_1d[k * P_1d + j];
1089       laplace[i + j * P_1d] = sum;
1090     }
1091   }
1092   CeedScalar perturbation = dim > 2 ? 1e-6 : 1e-4;
1093   for (CeedInt i = 0; i < P_1d; i++) laplace[i + P_1d * i] += perturbation;
1094   return CEED_ERROR_SUCCESS;
1095 }
1096 CeedPragmaOptimizeOn
1097 
1098 /// @}
1099 
1100 /// ----------------------------------------------------------------------------
1101 /// CeedOperator Backend API
1102 /// ----------------------------------------------------------------------------
1103 /// @addtogroup CeedOperatorBackend
1104 /// @{
1105 
1106 /**
1107   @brief Select correct basis matrix pointer based on @ref CeedEvalMode
1108 
1109   @param[in]  basis     `CeedBasis` from which to get the basis matrix
1110   @param[in]  eval_mode Current basis evaluation mode
1111   @param[in]  identity  Pointer to identity matrix
1112   @param[out] basis_ptr `CeedBasis` pointer to set
1113 
1114   @ref Backend
1115 **/
1116 int CeedOperatorGetBasisPointer(CeedBasis basis, CeedEvalMode eval_mode, const CeedScalar *identity, const CeedScalar **basis_ptr) {
1117   switch (eval_mode) {
1118     case CEED_EVAL_NONE:
1119       *basis_ptr = identity;
1120       break;
1121     case CEED_EVAL_INTERP:
1122       CeedCall(CeedBasisGetInterp(basis, basis_ptr));
1123       break;
1124     case CEED_EVAL_GRAD:
1125       CeedCall(CeedBasisGetGrad(basis, basis_ptr));
1126       break;
1127     case CEED_EVAL_DIV:
1128       CeedCall(CeedBasisGetDiv(basis, basis_ptr));
1129       break;
1130     case CEED_EVAL_CURL:
1131       CeedCall(CeedBasisGetCurl(basis, basis_ptr));
1132       break;
1133     case CEED_EVAL_WEIGHT:
1134       break;  // Caught by QF Assembly
1135   }
1136   assert(*basis_ptr != NULL);
1137   return CEED_ERROR_SUCCESS;
1138 }
1139 
1140 /**
1141   @brief Create point block restriction for active `CeedOperatorField`
1142 
1143   @param[in]  rstr             Original `CeedElemRestriction` for active field
1144   @param[out] point_block_rstr Address of the variable where the newly created `CeedElemRestriction` will be stored
1145 
1146   @return An error code: 0 - success, otherwise - failure
1147 
1148   @ref Backend
1149 **/
1150 int CeedOperatorCreateActivePointBlockRestriction(CeedElemRestriction rstr, CeedElemRestriction *point_block_rstr) {
1151   Ceed           ceed;
1152   CeedInt        num_elem, num_comp, shift, elem_size, comp_stride, *point_block_offsets;
1153   CeedSize       l_size;
1154   const CeedInt *offsets;
1155 
1156   CeedCall(CeedElemRestrictionGetCeed(rstr, &ceed));
1157   CeedCall(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets));
1158 
1159   // Expand offsets
1160   CeedCall(CeedElemRestrictionGetNumElements(rstr, &num_elem));
1161   CeedCall(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
1162   CeedCall(CeedElemRestrictionGetElementSize(rstr, &elem_size));
1163   CeedCall(CeedElemRestrictionGetCompStride(rstr, &comp_stride));
1164   CeedCall(CeedElemRestrictionGetLVectorSize(rstr, &l_size));
1165   shift = num_comp;
1166   if (comp_stride != 1) shift *= num_comp;
1167   CeedCall(CeedCalloc(num_elem * elem_size, &point_block_offsets));
1168   for (CeedInt i = 0; i < num_elem * elem_size; i++) {
1169     point_block_offsets[i] = offsets[i] * shift;
1170   }
1171 
1172   // Create new restriction
1173   CeedCall(CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp * num_comp, 1, l_size * num_comp, CEED_MEM_HOST, CEED_OWN_POINTER,
1174                                      point_block_offsets, point_block_rstr));
1175 
1176   // Cleanup
1177   CeedCall(CeedElemRestrictionRestoreOffsets(rstr, &offsets));
1178   CeedCall(CeedDestroy(&ceed));
1179   return CEED_ERROR_SUCCESS;
1180 }
1181 
1182 /**
1183   @brief Get `CeedQFunctionAssemblyData`
1184 
1185   @param[in]  op   `CeedOperator` to assemble
1186   @param[out] data `CeedQFunctionAssemblyData`
1187 
1188   @return An error code: 0 - success, otherwise - failure
1189 
1190   @ref Backend
1191 **/
1192 int CeedOperatorGetQFunctionAssemblyData(CeedOperator op, CeedQFunctionAssemblyData *data) {
1193   if (!op->qf_assembled) {
1194     CeedQFunctionAssemblyData data;
1195 
1196     CeedCall(CeedQFunctionAssemblyDataCreate(op->ceed, &data));
1197     op->qf_assembled = data;
1198   }
1199   *data = op->qf_assembled;
1200   return CEED_ERROR_SUCCESS;
1201 }
1202 
1203 /**
1204   @brief Create object holding `CeedQFunction` assembly data for `CeedOperator`
1205 
1206   @param[in]  ceed `Ceed` object used to create the `CeedQFunctionAssemblyData`
1207   @param[out] data Address of the variable where the newly created `CeedQFunctionAssemblyData` will be stored
1208 
1209   @return An error code: 0 - success, otherwise - failure
1210 
1211   @ref Backend
1212 **/
1213 int CeedQFunctionAssemblyDataCreate(Ceed ceed, CeedQFunctionAssemblyData *data) {
1214   CeedCall(CeedCalloc(1, data));
1215   (*data)->ref_count = 1;
1216   (*data)->ceed      = ceed;
1217   CeedCall(CeedReference(ceed));
1218   return CEED_ERROR_SUCCESS;
1219 }
1220 
1221 /**
1222   @brief Increment the reference counter for a `CeedQFunctionAssemblyData`
1223 
1224   @param[in,out] data `CeedQFunctionAssemblyData` to increment the reference counter
1225 
1226   @return An error code: 0 - success, otherwise - failure
1227 
1228   @ref Backend
1229 **/
1230 int CeedQFunctionAssemblyDataReference(CeedQFunctionAssemblyData data) {
1231   data->ref_count++;
1232   return CEED_ERROR_SUCCESS;
1233 }
1234 
1235 /**
1236   @brief Set re-use of `CeedQFunctionAssemblyData`
1237 
1238   @param[in,out] data       `CeedQFunctionAssemblyData` to mark for reuse
1239   @param[in]     reuse_data Boolean flag indicating data re-use
1240 
1241   @return An error code: 0 - success, otherwise - failure
1242 
1243   @ref Backend
1244 **/
1245 int CeedQFunctionAssemblyDataSetReuse(CeedQFunctionAssemblyData data, bool reuse_data) {
1246   data->reuse_data        = reuse_data;
1247   data->needs_data_update = true;
1248   return CEED_ERROR_SUCCESS;
1249 }
1250 
1251 /**
1252   @brief Mark `CeedQFunctionAssemblyData` as stale
1253 
1254   @param[in,out] data              `CeedQFunctionAssemblyData` to mark as stale
1255   @param[in]     needs_data_update Boolean flag indicating if update is needed or completed
1256 
1257   @return An error code: 0 - success, otherwise - failure
1258 
1259   @ref Backend
1260 **/
1261 int CeedQFunctionAssemblyDataSetUpdateNeeded(CeedQFunctionAssemblyData data, bool needs_data_update) {
1262   data->needs_data_update = needs_data_update;
1263   return CEED_ERROR_SUCCESS;
1264 }
1265 
1266 /**
1267   @brief Determine if `CeedQFunctionAssemblyData` needs update
1268 
1269   @param[in]  data             `CeedQFunctionAssemblyData` to mark as stale
1270   @param[out] is_update_needed Boolean flag indicating if re-assembly is required
1271 
1272   @return An error code: 0 - success, otherwise - failure
1273 
1274   @ref Backend
1275 **/
1276 int CeedQFunctionAssemblyDataIsUpdateNeeded(CeedQFunctionAssemblyData data, bool *is_update_needed) {
1277   *is_update_needed = !data->reuse_data || data->needs_data_update;
1278   return CEED_ERROR_SUCCESS;
1279 }
1280 
1281 /**
1282   @brief Copy the pointer to a `CeedQFunctionAssemblyData`.
1283 
1284   Both pointers should be destroyed with @ref CeedQFunctionAssemblyDataDestroy().
1285 
1286   Note: If the value of ` *data_copy` passed to this function is non-`NULL` , then it is assumed that ` *data_copy` is a pointer to a `CeedQFunctionAssemblyData`.
1287         This `CeedQFunctionAssemblyData` will be destroyed if ` *data_copy` is the only reference to this `CeedQFunctionAssemblyData`.
1288 
1289   @param[in]     data      `CeedQFunctionAssemblyData` to copy reference to
1290   @param[in,out] data_copy Variable to store copied reference
1291 
1292   @return An error code: 0 - success, otherwise - failure
1293 
1294   @ref Backend
1295 **/
1296 int CeedQFunctionAssemblyDataReferenceCopy(CeedQFunctionAssemblyData data, CeedQFunctionAssemblyData *data_copy) {
1297   CeedCall(CeedQFunctionAssemblyDataReference(data));
1298   CeedCall(CeedQFunctionAssemblyDataDestroy(data_copy));
1299   *data_copy = data;
1300   return CEED_ERROR_SUCCESS;
1301 }
1302 
1303 /**
1304   @brief Get setup status for internal objects for `CeedQFunctionAssemblyData`
1305 
1306   @param[in]  data     `CeedQFunctionAssemblyData` to retrieve status
1307   @param[out] is_setup Boolean flag for setup status
1308 
1309   @return An error code: 0 - success, otherwise - failure
1310 
1311   @ref Backend
1312 **/
1313 int CeedQFunctionAssemblyDataIsSetup(CeedQFunctionAssemblyData data, bool *is_setup) {
1314   *is_setup = data->is_setup;
1315   return CEED_ERROR_SUCCESS;
1316 }
1317 
1318 /**
1319   @brief Set internal objects for `CeedQFunctionAssemblyData`
1320 
1321   @param[in,out] data `CeedQFunctionAssemblyData` to set objects
1322   @param[in]     vec  `CeedVector` to store assembled `CeedQFunction` at quadrature points
1323   @param[in]     rstr `CeedElemRestriction` for `CeedVector` containing assembled `CeedQFunction`
1324 
1325   @return An error code: 0 - success, otherwise - failure
1326 
1327   @ref Backend
1328 **/
1329 int CeedQFunctionAssemblyDataSetObjects(CeedQFunctionAssemblyData data, CeedVector vec, CeedElemRestriction rstr) {
1330   CeedCall(CeedVectorReferenceCopy(vec, &data->vec));
1331   CeedCall(CeedElemRestrictionReferenceCopy(rstr, &data->rstr));
1332 
1333   data->is_setup = true;
1334   return CEED_ERROR_SUCCESS;
1335 }
1336 
1337 /**
1338   @brief Get internal objects for `CeedQFunctionAssemblyData`
1339 
1340   @param[in,out] data `CeedQFunctionAssemblyData` to set objects
1341   @param[out]    vec  `CeedVector` to store assembled `CeedQFunction` at quadrature points
1342   @param[out]    rstr `CeedElemRestriction` for `CeedVector` containing assembled `CeedQFunction`
1343 
1344   @return An error code: 0 - success, otherwise - failure
1345 
1346   @ref Backend
1347 **/
1348 int CeedQFunctionAssemblyDataGetObjects(CeedQFunctionAssemblyData data, CeedVector *vec, CeedElemRestriction *rstr) {
1349   CeedCheck(data->is_setup, data->ceed, CEED_ERROR_INCOMPLETE, "Internal objects not set; must call CeedQFunctionAssemblyDataSetObjects first.");
1350 
1351   CeedCall(CeedVectorReferenceCopy(data->vec, vec));
1352   CeedCall(CeedElemRestrictionReferenceCopy(data->rstr, rstr));
1353   return CEED_ERROR_SUCCESS;
1354 }
1355 
1356 /**
1357   @brief Destroy `CeedQFunctionAssemblyData`
1358 
1359   @param[in,out] data  `CeedQFunctionAssemblyData` to destroy
1360 
1361   @return An error code: 0 - success, otherwise - failure
1362 
1363   @ref Backend
1364 **/
1365 int CeedQFunctionAssemblyDataDestroy(CeedQFunctionAssemblyData *data) {
1366   if (!*data || --(*data)->ref_count > 0) {
1367     *data = NULL;
1368     return CEED_ERROR_SUCCESS;
1369   }
1370   CeedCall(CeedDestroy(&(*data)->ceed));
1371   CeedCall(CeedVectorDestroy(&(*data)->vec));
1372   CeedCall(CeedElemRestrictionDestroy(&(*data)->rstr));
1373 
1374   CeedCall(CeedFree(data));
1375   return CEED_ERROR_SUCCESS;
1376 }
1377 
1378 /**
1379   @brief Get `CeedOperatorAssemblyData`
1380 
1381   @param[in]  op   `CeedOperator` to assemble
1382   @param[out] data `CeedOperatorAssemblyData`
1383 
1384   @return An error code: 0 - success, otherwise - failure
1385 
1386   @ref Backend
1387 **/
1388 int CeedOperatorGetOperatorAssemblyData(CeedOperator op, CeedOperatorAssemblyData *data) {
1389   if (!op->op_assembled) {
1390     CeedOperatorAssemblyData data;
1391 
1392     CeedCall(CeedOperatorAssemblyDataCreate(op->ceed, op, &data));
1393     op->op_assembled = data;
1394   }
1395   *data = op->op_assembled;
1396   return CEED_ERROR_SUCCESS;
1397 }
1398 
1399 /**
1400   @brief Create object holding `CeedOperator` assembly data.
1401 
1402   The `CeedOperatorAssemblyData` holds an array with references to every active `CeedBasis` used in the `CeedOperator`.
1403   An array with references to the corresponding active `CeedElemRestriction` is also stored.
1404   For each active `CeedBasis, the `CeedOperatorAssemblyData` holds an array of all input and output @ref CeedEvalMode for this `CeedBasis`.
1405   The `CeedOperatorAssemblyData` holds an array of offsets for indexing into the assembled `CeedQFunction` arrays to the row representing each @ref CeedEvalMode.
1406   The number of input columns across all active bases for the assembled `CeedQFunction` is also stored.
1407   Lastly, the `CeedOperatorAssembly` data holds assembled matrices representing the full action of the `CeedBasis` for all @ref CeedEvalMode.
1408 
1409   @param[in]  ceed `Ceed` object used to create the `CeedOperatorAssemblyData`
1410   @param[in]  op   `CeedOperator` to be assembled
1411   @param[out] data Address of the variable where the newly created `CeedOperatorAssemblyData` will be stored
1412 
1413   @return An error code: 0 - success, otherwise - failure
1414 
1415   @ref Backend
1416 **/
1417 int CeedOperatorAssemblyDataCreate(Ceed ceed, CeedOperator op, CeedOperatorAssemblyData *data) {
1418   CeedInt             num_active_bases_in = 0, num_active_bases_out = 0, offset = 0;
1419   CeedInt             num_input_fields, *num_eval_modes_in = NULL, num_output_fields, *num_eval_modes_out = NULL;
1420   CeedSize          **eval_mode_offsets_in = NULL, **eval_mode_offsets_out = NULL;
1421   CeedEvalMode      **eval_modes_in = NULL, **eval_modes_out = NULL;
1422   CeedQFunctionField *qf_fields;
1423   CeedQFunction       qf;
1424   CeedOperatorField  *op_fields;
1425   bool                is_composite;
1426 
1427   CeedCall(CeedOperatorIsComposite(op, &is_composite));
1428   CeedCheck(!is_composite, ceed, CEED_ERROR_INCOMPATIBLE, "Can only create CeedOperator assembly data for non-composite operators.");
1429 
1430   // Allocate
1431   CeedCall(CeedCalloc(1, data));
1432   (*data)->ceed = ceed;
1433   CeedCall(CeedReference(ceed));
1434 
1435   // Build OperatorAssembly data
1436   CeedCall(CeedOperatorGetQFunction(op, &qf));
1437 
1438   // Determine active input basis
1439   CeedCall(CeedQFunctionGetFields(qf, &num_input_fields, &qf_fields, NULL, NULL));
1440   CeedCall(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1441   for (CeedInt i = 0; i < num_input_fields; i++) {
1442     CeedVector vec;
1443 
1444     CeedCall(CeedOperatorFieldGetVector(op_fields[i], &vec));
1445     if (vec == CEED_VECTOR_ACTIVE) {
1446       CeedInt      index = -1, num_comp, q_comp;
1447       CeedEvalMode eval_mode;
1448       CeedBasis    basis_in = NULL;
1449 
1450       CeedCall(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
1451       CeedCall(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1452       CeedCall(CeedBasisGetNumComponents(basis_in, &num_comp));
1453       CeedCall(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1454       for (CeedInt i = 0; i < num_active_bases_in; i++) {
1455         if ((*data)->active_bases_in[i] == basis_in) index = i;
1456       }
1457       if (index == -1) {
1458         CeedElemRestriction elem_rstr_in;
1459 
1460         index = num_active_bases_in;
1461         CeedCall(CeedRealloc(num_active_bases_in + 1, &(*data)->active_bases_in));
1462         (*data)->active_bases_in[num_active_bases_in] = NULL;
1463         CeedCall(CeedBasisReferenceCopy(basis_in, &(*data)->active_bases_in[num_active_bases_in]));
1464         CeedCall(CeedRealloc(num_active_bases_in + 1, &(*data)->active_elem_rstrs_in));
1465         (*data)->active_elem_rstrs_in[num_active_bases_in] = NULL;
1466         CeedCall(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr_in));
1467         CeedCall(CeedElemRestrictionReferenceCopy(elem_rstr_in, &(*data)->active_elem_rstrs_in[num_active_bases_in]));
1468         CeedCall(CeedElemRestrictionDestroy(&elem_rstr_in));
1469         CeedCall(CeedRealloc(num_active_bases_in + 1, &num_eval_modes_in));
1470         num_eval_modes_in[index] = 0;
1471         CeedCall(CeedRealloc(num_active_bases_in + 1, &eval_modes_in));
1472         eval_modes_in[index] = NULL;
1473         CeedCall(CeedRealloc(num_active_bases_in + 1, &eval_mode_offsets_in));
1474         eval_mode_offsets_in[index] = NULL;
1475         CeedCall(CeedRealloc(num_active_bases_in + 1, &(*data)->assembled_bases_in));
1476         (*data)->assembled_bases_in[index] = NULL;
1477         num_active_bases_in++;
1478       }
1479       if (eval_mode != CEED_EVAL_WEIGHT) {
1480         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1481         CeedCall(CeedRealloc(num_eval_modes_in[index] + q_comp, &eval_modes_in[index]));
1482         CeedCall(CeedRealloc(num_eval_modes_in[index] + q_comp, &eval_mode_offsets_in[index]));
1483         for (CeedInt d = 0; d < q_comp; d++) {
1484           eval_modes_in[index][num_eval_modes_in[index] + d]        = eval_mode;
1485           eval_mode_offsets_in[index][num_eval_modes_in[index] + d] = offset;
1486           offset += num_comp;
1487         }
1488         num_eval_modes_in[index] += q_comp;
1489       }
1490       CeedCall(CeedBasisDestroy(&basis_in));
1491     }
1492     CeedCall(CeedVectorDestroy(&vec));
1493   }
1494 
1495   // Determine active output basis
1496   CeedCall(CeedQFunctionGetFields(qf, NULL, NULL, &num_output_fields, &qf_fields));
1497   CeedCall(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1498   offset = 0;
1499   for (CeedInt i = 0; i < num_output_fields; i++) {
1500     CeedVector vec;
1501 
1502     CeedCall(CeedOperatorFieldGetVector(op_fields[i], &vec));
1503     if (vec == CEED_VECTOR_ACTIVE) {
1504       CeedInt      index = -1, num_comp, q_comp;
1505       CeedEvalMode eval_mode;
1506       CeedBasis    basis_out = NULL;
1507 
1508       CeedCall(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
1509       CeedCall(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1510       CeedCall(CeedBasisGetNumComponents(basis_out, &num_comp));
1511       CeedCall(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1512       for (CeedInt i = 0; i < num_active_bases_out; i++) {
1513         if ((*data)->active_bases_out[i] == basis_out) index = i;
1514       }
1515       if (index == -1) {
1516         CeedElemRestriction elem_rstr_out;
1517 
1518         index = num_active_bases_out;
1519         CeedCall(CeedRealloc(num_active_bases_out + 1, &(*data)->active_bases_out));
1520         (*data)->active_bases_out[num_active_bases_out] = NULL;
1521         CeedCall(CeedBasisReferenceCopy(basis_out, &(*data)->active_bases_out[num_active_bases_out]));
1522         CeedCall(CeedRealloc(num_active_bases_out + 1, &(*data)->active_elem_rstrs_out));
1523         (*data)->active_elem_rstrs_out[num_active_bases_out] = NULL;
1524         CeedCall(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr_out));
1525         CeedCall(CeedElemRestrictionReferenceCopy(elem_rstr_out, &(*data)->active_elem_rstrs_out[num_active_bases_out]));
1526         CeedCall(CeedElemRestrictionDestroy(&elem_rstr_out));
1527         CeedCall(CeedRealloc(num_active_bases_out + 1, &num_eval_modes_out));
1528         num_eval_modes_out[index] = 0;
1529         CeedCall(CeedRealloc(num_active_bases_out + 1, &eval_modes_out));
1530         eval_modes_out[index] = NULL;
1531         CeedCall(CeedRealloc(num_active_bases_out + 1, &eval_mode_offsets_out));
1532         eval_mode_offsets_out[index] = NULL;
1533         CeedCall(CeedRealloc(num_active_bases_out + 1, &(*data)->assembled_bases_out));
1534         (*data)->assembled_bases_out[index] = NULL;
1535         num_active_bases_out++;
1536       }
1537       if (eval_mode != CEED_EVAL_WEIGHT) {
1538         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1539         CeedCall(CeedRealloc(num_eval_modes_out[index] + q_comp, &eval_modes_out[index]));
1540         CeedCall(CeedRealloc(num_eval_modes_out[index] + q_comp, &eval_mode_offsets_out[index]));
1541         for (CeedInt d = 0; d < q_comp; d++) {
1542           eval_modes_out[index][num_eval_modes_out[index] + d]        = eval_mode;
1543           eval_mode_offsets_out[index][num_eval_modes_out[index] + d] = offset;
1544           offset += num_comp;
1545         }
1546         num_eval_modes_out[index] += q_comp;
1547       }
1548       CeedCall(CeedBasisDestroy(&basis_out));
1549     }
1550     CeedCall(CeedVectorDestroy(&vec));
1551   }
1552   CeedCall(CeedQFunctionDestroy(&qf));
1553   (*data)->num_active_bases_in   = num_active_bases_in;
1554   (*data)->num_eval_modes_in     = num_eval_modes_in;
1555   (*data)->eval_modes_in         = eval_modes_in;
1556   (*data)->eval_mode_offsets_in  = eval_mode_offsets_in;
1557   (*data)->num_active_bases_out  = num_active_bases_out;
1558   (*data)->num_eval_modes_out    = num_eval_modes_out;
1559   (*data)->eval_modes_out        = eval_modes_out;
1560   (*data)->eval_mode_offsets_out = eval_mode_offsets_out;
1561   (*data)->num_output_components = offset;
1562   return CEED_ERROR_SUCCESS;
1563 }
1564 
1565 /**
1566   @brief Get `CeedOperator` @ref CeedEvalMode for assembly.
1567 
1568   Note: See @ref CeedOperatorAssemblyDataCreate() for a full description of the data stored in this object.
1569 
1570   @param[in]  data                  `CeedOperatorAssemblyData`
1571   @param[out] num_active_bases_in   Total number of active bases for input
1572   @param[out] num_eval_modes_in     Pointer to hold array of numbers of input @ref CeedEvalMode, or `NULL`.
1573                                       `eval_modes_in[0]` holds an array of eval modes for the first active `CeedBasis`.
1574   @param[out] eval_modes_in         Pointer to hold arrays of input @ref CeedEvalMode, or `NULL`
1575   @param[out] eval_mode_offsets_in  Pointer to hold arrays of input offsets at each quadrature point
1576   @param[out] num_active_bases_out  Total number of active bases for output
1577   @param[out] num_eval_modes_out    Pointer to hold array of numbers of output @ref CeedEvalMode, or `NULL`
1578   @param[out] eval_modes_out        Pointer to hold arrays of output @ref CeedEvalMode, or `NULL`
1579   @param[out] eval_mode_offsets_out Pointer to hold arrays of output offsets at each quadrature point
1580   @param[out] num_output_components The number of columns in the assembled `CeedQFunction` matrix for each quadrature point, including contributions of all active bases
1581 
1582   @return An error code: 0 - success, otherwise - failure
1583 
1584   @ref Backend
1585 **/
1586 int CeedOperatorAssemblyDataGetEvalModes(CeedOperatorAssemblyData data, CeedInt *num_active_bases_in, CeedInt **num_eval_modes_in,
1587                                          const CeedEvalMode ***eval_modes_in, CeedSize ***eval_mode_offsets_in, CeedInt *num_active_bases_out,
1588                                          CeedInt **num_eval_modes_out, const CeedEvalMode ***eval_modes_out, CeedSize ***eval_mode_offsets_out,
1589                                          CeedSize *num_output_components) {
1590   if (num_active_bases_in) *num_active_bases_in = data->num_active_bases_in;
1591   if (num_eval_modes_in) *num_eval_modes_in = data->num_eval_modes_in;
1592   if (eval_modes_in) *eval_modes_in = (const CeedEvalMode **)data->eval_modes_in;
1593   if (eval_mode_offsets_in) *eval_mode_offsets_in = data->eval_mode_offsets_in;
1594   if (num_active_bases_out) *num_active_bases_out = data->num_active_bases_out;
1595   if (num_eval_modes_out) *num_eval_modes_out = data->num_eval_modes_out;
1596   if (eval_modes_out) *eval_modes_out = (const CeedEvalMode **)data->eval_modes_out;
1597   if (eval_mode_offsets_out) *eval_mode_offsets_out = data->eval_mode_offsets_out;
1598   if (num_output_components) *num_output_components = data->num_output_components;
1599   return CEED_ERROR_SUCCESS;
1600 }
1601 
1602 /**
1603   @brief Get `CeedOperator` `CeedBasis` data for assembly.
1604 
1605   Note: See @ref CeedOperatorAssemblyDataCreate() for a full description of the data stored in this object.
1606 
1607   @param[in]  data                 `CeedOperatorAssemblyData`
1608   @param[out] num_active_bases_in  Number of active input bases, or `NULL`
1609   @param[out] active_bases_in      Pointer to hold active input `CeedBasis`, or `NULL`
1610   @param[out] assembled_bases_in   Pointer to hold assembled active input `B` , or `NULL`
1611   @param[out] num_active_bases_out Number of active output bases, or `NULL`
1612   @param[out] active_bases_out     Pointer to hold active output `CeedBasis`, or `NULL`
1613   @param[out] assembled_bases_out  Pointer to hold assembled active output `B` , or `NULL`
1614 
1615   @return An error code: 0 - success, otherwise - failure
1616 
1617   @ref Backend
1618 **/
1619 int CeedOperatorAssemblyDataGetBases(CeedOperatorAssemblyData data, CeedInt *num_active_bases_in, CeedBasis **active_bases_in,
1620                                      const CeedScalar ***assembled_bases_in, CeedInt *num_active_bases_out, CeedBasis **active_bases_out,
1621                                      const CeedScalar ***assembled_bases_out) {
1622   // Assemble B_in, B_out if needed
1623   if (assembled_bases_in && !data->assembled_bases_in[0]) {
1624     CeedInt num_qpts;
1625 
1626     if (data->active_bases_in[0] == CEED_BASIS_NONE) CeedCall(CeedElemRestrictionGetElementSize(data->active_elem_rstrs_in[0], &num_qpts));
1627     else CeedCall(CeedBasisGetNumQuadraturePoints(data->active_bases_in[0], &num_qpts));
1628     for (CeedInt b = 0; b < data->num_active_bases_in; b++) {
1629       bool        has_eval_none = false;
1630       CeedInt     num_nodes;
1631       CeedScalar *B_in = NULL, *identity = NULL;
1632 
1633       CeedCall(CeedElemRestrictionGetElementSize(data->active_elem_rstrs_in[b], &num_nodes));
1634       CeedCall(CeedCalloc(num_qpts * num_nodes * data->num_eval_modes_in[b], &B_in));
1635 
1636       for (CeedInt i = 0; i < data->num_eval_modes_in[b]; i++) {
1637         has_eval_none = has_eval_none || (data->eval_modes_in[b][i] == CEED_EVAL_NONE);
1638       }
1639       if (has_eval_none) {
1640         CeedCall(CeedCalloc(num_qpts * num_nodes, &identity));
1641         for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) {
1642           identity[i * num_nodes + i] = 1.0;
1643         }
1644       }
1645 
1646       for (CeedInt q = 0; q < num_qpts; q++) {
1647         for (CeedInt n = 0; n < num_nodes; n++) {
1648           CeedInt      d_in              = 0, q_comp_in;
1649           CeedEvalMode eval_mode_in_prev = CEED_EVAL_NONE;
1650 
1651           for (CeedInt e_in = 0; e_in < data->num_eval_modes_in[b]; e_in++) {
1652             const CeedInt     qq = data->num_eval_modes_in[b] * q;
1653             const CeedScalar *B  = NULL;
1654 
1655             CeedCall(CeedOperatorGetBasisPointer(data->active_bases_in[b], data->eval_modes_in[b][e_in], identity, &B));
1656             CeedCall(CeedBasisGetNumQuadratureComponents(data->active_bases_in[b], data->eval_modes_in[b][e_in], &q_comp_in));
1657             if (q_comp_in > 1) {
1658               if (e_in == 0 || data->eval_modes_in[b][e_in] != eval_mode_in_prev) d_in = 0;
1659               else B = &B[(++d_in) * num_qpts * num_nodes];
1660             }
1661             eval_mode_in_prev                 = data->eval_modes_in[b][e_in];
1662             B_in[(qq + e_in) * num_nodes + n] = B[q * num_nodes + n];
1663           }
1664         }
1665       }
1666       if (identity) CeedCall(CeedFree(&identity));
1667       data->assembled_bases_in[b] = B_in;
1668     }
1669   }
1670 
1671   if (assembled_bases_out && !data->assembled_bases_out[0]) {
1672     CeedInt num_qpts;
1673 
1674     if (data->active_bases_out[0] == CEED_BASIS_NONE) CeedCall(CeedElemRestrictionGetElementSize(data->active_elem_rstrs_out[0], &num_qpts));
1675     else CeedCall(CeedBasisGetNumQuadraturePoints(data->active_bases_out[0], &num_qpts));
1676     for (CeedInt b = 0; b < data->num_active_bases_out; b++) {
1677       bool        has_eval_none = false;
1678       CeedInt     num_nodes;
1679       CeedScalar *B_out = NULL, *identity = NULL;
1680 
1681       CeedCall(CeedElemRestrictionGetElementSize(data->active_elem_rstrs_out[b], &num_nodes));
1682       CeedCall(CeedCalloc(num_qpts * num_nodes * data->num_eval_modes_out[b], &B_out));
1683 
1684       for (CeedInt i = 0; i < data->num_eval_modes_out[b]; i++) {
1685         has_eval_none = has_eval_none || (data->eval_modes_out[b][i] == CEED_EVAL_NONE);
1686       }
1687       if (has_eval_none) {
1688         CeedCall(CeedCalloc(num_qpts * num_nodes, &identity));
1689         for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) {
1690           identity[i * num_nodes + i] = 1.0;
1691         }
1692       }
1693 
1694       for (CeedInt q = 0; q < num_qpts; q++) {
1695         for (CeedInt n = 0; n < num_nodes; n++) {
1696           CeedInt      d_out              = 0, q_comp_out;
1697           CeedEvalMode eval_mode_out_prev = CEED_EVAL_NONE;
1698 
1699           for (CeedInt e_out = 0; e_out < data->num_eval_modes_out[b]; e_out++) {
1700             const CeedInt     qq = data->num_eval_modes_out[b] * q;
1701             const CeedScalar *B  = NULL;
1702 
1703             CeedCall(CeedOperatorGetBasisPointer(data->active_bases_out[b], data->eval_modes_out[b][e_out], identity, &B));
1704             CeedCall(CeedBasisGetNumQuadratureComponents(data->active_bases_out[b], data->eval_modes_out[b][e_out], &q_comp_out));
1705             if (q_comp_out > 1) {
1706               if (e_out == 0 || data->eval_modes_out[b][e_out] != eval_mode_out_prev) d_out = 0;
1707               else B = &B[(++d_out) * num_qpts * num_nodes];
1708             }
1709             eval_mode_out_prev                  = data->eval_modes_out[b][e_out];
1710             B_out[(qq + e_out) * num_nodes + n] = B[q * num_nodes + n];
1711           }
1712         }
1713       }
1714       if (identity) CeedCall(CeedFree(&identity));
1715       data->assembled_bases_out[b] = B_out;
1716     }
1717   }
1718 
1719   // Pass out assembled data
1720   if (num_active_bases_in) *num_active_bases_in = data->num_active_bases_in;
1721   if (active_bases_in) *active_bases_in = data->active_bases_in;
1722   if (assembled_bases_in) *assembled_bases_in = (const CeedScalar **)data->assembled_bases_in;
1723   if (num_active_bases_out) *num_active_bases_out = data->num_active_bases_out;
1724   if (active_bases_out) *active_bases_out = data->active_bases_out;
1725   if (assembled_bases_out) *assembled_bases_out = (const CeedScalar **)data->assembled_bases_out;
1726   return CEED_ERROR_SUCCESS;
1727 }
1728 
1729 /**
1730   @brief Get `CeedOperator` `CeedBasis` data for assembly.
1731 
1732   Note: See @ref CeedOperatorAssemblyDataCreate() for a full description of the data stored in this object.
1733 
1734   @param[in]  data                      `CeedOperatorAssemblyData`
1735   @param[out] num_active_elem_rstrs_in  Number of active input element restrictions, or `NULL`
1736   @param[out] active_elem_rstrs_in      Pointer to hold active input `CeedElemRestriction`, or `NULL`
1737   @param[out] num_active_elem_rstrs_out Number of active output element restrictions, or `NULL`
1738   @param[out] active_elem_rstrs_out     Pointer to hold active output `CeedElemRestriction`, or `NULL`
1739 
1740   @return An error code: 0 - success, otherwise - failure
1741 
1742   @ref Backend
1743 **/
1744 int CeedOperatorAssemblyDataGetElemRestrictions(CeedOperatorAssemblyData data, CeedInt *num_active_elem_rstrs_in,
1745                                                 CeedElemRestriction **active_elem_rstrs_in, CeedInt *num_active_elem_rstrs_out,
1746                                                 CeedElemRestriction **active_elem_rstrs_out) {
1747   if (num_active_elem_rstrs_in) *num_active_elem_rstrs_in = data->num_active_bases_in;
1748   if (active_elem_rstrs_in) *active_elem_rstrs_in = data->active_elem_rstrs_in;
1749   if (num_active_elem_rstrs_out) *num_active_elem_rstrs_out = data->num_active_bases_out;
1750   if (active_elem_rstrs_out) *active_elem_rstrs_out = data->active_elem_rstrs_out;
1751   return CEED_ERROR_SUCCESS;
1752 }
1753 
1754 /**
1755   @brief Destroy `CeedOperatorAssemblyData`
1756 
1757   @param[in,out] data `CeedOperatorAssemblyData` to destroy
1758 
1759   @return An error code: 0 - success, otherwise - failure
1760 
1761   @ref Backend
1762 **/
1763 int CeedOperatorAssemblyDataDestroy(CeedOperatorAssemblyData *data) {
1764   if (!*data) {
1765     *data = NULL;
1766     return CEED_ERROR_SUCCESS;
1767   }
1768   CeedCall(CeedDestroy(&(*data)->ceed));
1769   for (CeedInt b = 0; b < (*data)->num_active_bases_in; b++) {
1770     CeedCall(CeedBasisDestroy(&(*data)->active_bases_in[b]));
1771     CeedCall(CeedElemRestrictionDestroy(&(*data)->active_elem_rstrs_in[b]));
1772     CeedCall(CeedFree(&(*data)->eval_modes_in[b]));
1773     CeedCall(CeedFree(&(*data)->eval_mode_offsets_in[b]));
1774     CeedCall(CeedFree(&(*data)->assembled_bases_in[b]));
1775   }
1776   for (CeedInt b = 0; b < (*data)->num_active_bases_out; b++) {
1777     CeedCall(CeedBasisDestroy(&(*data)->active_bases_out[b]));
1778     CeedCall(CeedElemRestrictionDestroy(&(*data)->active_elem_rstrs_out[b]));
1779     CeedCall(CeedFree(&(*data)->eval_modes_out[b]));
1780     CeedCall(CeedFree(&(*data)->eval_mode_offsets_out[b]));
1781     CeedCall(CeedFree(&(*data)->assembled_bases_out[b]));
1782   }
1783   CeedCall(CeedFree(&(*data)->active_bases_in));
1784   CeedCall(CeedFree(&(*data)->active_bases_out));
1785   CeedCall(CeedFree(&(*data)->active_elem_rstrs_in));
1786   CeedCall(CeedFree(&(*data)->active_elem_rstrs_out));
1787   CeedCall(CeedFree(&(*data)->num_eval_modes_in));
1788   CeedCall(CeedFree(&(*data)->num_eval_modes_out));
1789   CeedCall(CeedFree(&(*data)->eval_modes_in));
1790   CeedCall(CeedFree(&(*data)->eval_modes_out));
1791   CeedCall(CeedFree(&(*data)->eval_mode_offsets_in));
1792   CeedCall(CeedFree(&(*data)->eval_mode_offsets_out));
1793   CeedCall(CeedFree(&(*data)->assembled_bases_in));
1794   CeedCall(CeedFree(&(*data)->assembled_bases_out));
1795 
1796   CeedCall(CeedFree(data));
1797   return CEED_ERROR_SUCCESS;
1798 }
1799 
1800 /**
1801   @brief Retrieve fallback `CeedOperator` with a reference `Ceed` for advanced `CeedOperator` functionality
1802 
1803   @param[in]  op          `CeedOperator` to retrieve fallback for
1804   @param[out] op_fallback Fallback `CeedOperator`
1805 
1806   @return An error code: 0 - success, otherwise - failure
1807 
1808   @ref Backend
1809 **/
1810 int CeedOperatorGetFallback(CeedOperator op, CeedOperator *op_fallback) {
1811   // Create if needed
1812   if (!op->op_fallback) CeedCall(CeedOperatorCreateFallback(op));
1813   if (op->op_fallback) {
1814     bool is_debug;
1815     Ceed ceed;
1816 
1817     CeedCall(CeedOperatorGetCeed(op, &ceed));
1818     CeedCall(CeedIsDebug(ceed, &is_debug));
1819     if (is_debug) {
1820       Ceed        ceed_fallback;
1821       const char *resource, *resource_fallback;
1822 
1823       CeedCall(CeedGetOperatorFallbackCeed(ceed, &ceed_fallback));
1824       CeedCall(CeedGetResource(ceed, &resource));
1825       CeedCall(CeedGetResource(ceed_fallback, &resource_fallback));
1826 
1827       CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "---------- CeedOperator Fallback ----------\n");
1828       CeedDebug(ceed, "Falling back from %s operator at address %p to %s operator at address %p\n", resource, op, resource_fallback, op->op_fallback);
1829       CeedCall(CeedDestroy(&ceed_fallback));
1830     }
1831     CeedCall(CeedDestroy(&ceed));
1832   }
1833   *op_fallback = op->op_fallback;
1834   return CEED_ERROR_SUCCESS;
1835 }
1836 
1837 /**
1838   @brief Get the parent `CeedOperator` for a fallback `CeedOperator`
1839 
1840   @param[in]  op     `CeedOperator` context
1841   @param[out] parent Variable to store parent `CeedOperator` context
1842 
1843   @return An error code: 0 - success, otherwise - failure
1844 
1845   @ref Backend
1846 **/
1847 int CeedOperatorGetFallbackParent(CeedOperator op, CeedOperator *parent) {
1848   *parent = op->op_fallback_parent ? op->op_fallback_parent : NULL;
1849   return CEED_ERROR_SUCCESS;
1850 }
1851 
1852 /**
1853   @brief Get the `Ceed` context of the parent `CeedOperator` for a fallback `CeedOperator`
1854 
1855   @param[in]  op     `CeedOperator` context
1856   @param[out] parent Variable to store parent `Ceed` context
1857 
1858   @return An error code: 0 - success, otherwise - failure
1859 
1860   @ref Backend
1861 **/
1862 int CeedOperatorGetFallbackParentCeed(CeedOperator op, Ceed *parent) {
1863   *parent = NULL;
1864   if (op->op_fallback_parent) CeedCall(CeedReferenceCopy(op->op_fallback_parent->ceed, parent));
1865   else CeedCall(CeedReferenceCopy(CeedOperatorReturnCeed(op), parent));
1866   return CEED_ERROR_SUCCESS;
1867 }
1868 
1869 /// @}
1870 
1871 /// ----------------------------------------------------------------------------
1872 /// CeedOperator Public API
1873 /// ----------------------------------------------------------------------------
1874 /// @addtogroup CeedOperatorUser
1875 /// @{
1876 
1877 /**
1878   @brief Assemble a linear `CeedQFunction` associated with a `CeedOperator`.
1879 
1880   This returns a `CeedVector` containing a matrix at each quadrature point providing the action of the `CeedQFunction` associated with the `CeedOperator`.
1881   The vector `assembled` is of shape `[num_elements, num_input_fields, num_output_fields, num_quad_points]` and contains column-major matrices representing the action of the `CeedQFunction` for a corresponding quadrature point on an element.
1882 
1883   Inputs and outputs are in the order provided by the user when adding `CeedOperator` fields.
1884   For example, a `CeedQFunction` with inputs `u` and `gradu` and outputs `gradv` and `v` , provided in that order, would result in an assembled `CeedQFunction` that consists of `(1 + dim) x (dim + 1)` matrices at each quadrature point acting on the input ` [u, du_0, du_1]` and producing the output `[dv_0, dv_1, v]`.
1885 
1886   Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
1887 
1888   @param[in]  op        `CeedOperator` to assemble `CeedQFunction`
1889   @param[out] assembled `CeedVector` to store assembled `CeedQFunction` at quadrature points
1890   @param[out] rstr      `CeedElemRestriction` for `CeedVector` containing assembled `CeedQFunction`
1891   @param[in]  request   Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
1892 
1893   @return An error code: 0 - success, otherwise - failure
1894 
1895   @ref User
1896 **/
1897 int CeedOperatorLinearAssembleQFunction(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1898   CeedCall(CeedOperatorCheckReady(op));
1899 
1900   if (op->LinearAssembleQFunction) {
1901     // Backend version
1902     CeedCall(op->LinearAssembleQFunction(op, assembled, rstr, request));
1903   } else {
1904     // Operator fallback
1905     CeedOperator op_fallback;
1906 
1907     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
1908     if (op_fallback) CeedCall(CeedOperatorLinearAssembleQFunction(op_fallback, assembled, rstr, request));
1909     else return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedOperatorLinearAssembleQFunction");
1910   }
1911   return CEED_ERROR_SUCCESS;
1912 }
1913 
1914 /**
1915   @brief Assemble `CeedQFunction` and store result internally.
1916 
1917   Return copied references of stored data to the caller.
1918   Caller is responsible for ownership and destruction of the copied references.
1919   See also @ref CeedOperatorLinearAssembleQFunction().
1920 
1921   Note: If the value of `assembled` or `rstr` passed to this function are non-`NULL` , then it is assumed that they hold valid pointers.
1922         These objects will be destroyed if `*assembled` or `*rstr` is the only reference to the object.
1923 
1924   @param[in]  op        `CeedOperator` to assemble `CeedQFunction`
1925   @param[out] assembled `CeedVector` to store assembled `CeedQFunction` at quadrature points
1926   @param[out] rstr      `CeedElemRestriction` for `CeedVector` containing assembled `CeedQFunction`
1927   @param[in]  request   Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
1928 
1929   @return An error code: 0 - success, otherwise - failure
1930 
1931   @ref User
1932 **/
1933 int CeedOperatorLinearAssembleQFunctionBuildOrUpdate(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1934   int (*LinearAssembleQFunctionUpdate)(CeedOperator, CeedVector, CeedElemRestriction, CeedRequest *) = NULL;
1935   CeedOperator op_assemble                                                                           = NULL;
1936   CeedOperator op_fallback_parent                                                                    = NULL;
1937 
1938   CeedCall(CeedOperatorCheckReady(op));
1939 
1940   // Determine if fallback parent or operator has implementation
1941   CeedCall(CeedOperatorGetFallbackParent(op, &op_fallback_parent));
1942   if (op_fallback_parent && op_fallback_parent->LinearAssembleQFunctionUpdate) {
1943     // -- Backend version for op fallback parent is faster, if it exists
1944     LinearAssembleQFunctionUpdate = op_fallback_parent->LinearAssembleQFunctionUpdate;
1945     op_assemble                   = op_fallback_parent;
1946   } else if (op->LinearAssembleQFunctionUpdate) {
1947     // -- Backend version for op
1948     LinearAssembleQFunctionUpdate = op->LinearAssembleQFunctionUpdate;
1949     op_assemble                   = op;
1950   }
1951 
1952   // Assemble QFunction
1953   if (LinearAssembleQFunctionUpdate) {
1954     // Backend or fallback parent version
1955     CeedQFunctionAssemblyData data;
1956     bool                      data_is_setup;
1957     CeedVector                assembled_vec  = NULL;
1958     CeedElemRestriction       assembled_rstr = NULL;
1959 
1960     CeedCall(CeedOperatorGetQFunctionAssemblyData(op, &data));
1961     CeedCall(CeedQFunctionAssemblyDataIsSetup(data, &data_is_setup));
1962     if (data_is_setup) {
1963       bool update_needed;
1964 
1965       CeedCall(CeedQFunctionAssemblyDataGetObjects(data, &assembled_vec, &assembled_rstr));
1966       CeedCall(CeedQFunctionAssemblyDataIsUpdateNeeded(data, &update_needed));
1967       if (update_needed) CeedCall(LinearAssembleQFunctionUpdate(op_assemble, assembled_vec, assembled_rstr, request));
1968     } else {
1969       CeedCall(CeedOperatorLinearAssembleQFunction(op_assemble, &assembled_vec, &assembled_rstr, request));
1970       CeedCall(CeedQFunctionAssemblyDataSetObjects(data, assembled_vec, assembled_rstr));
1971     }
1972     CeedCall(CeedQFunctionAssemblyDataSetUpdateNeeded(data, false));
1973 
1974     // Copy reference from internally held copy
1975     CeedCall(CeedVectorReferenceCopy(assembled_vec, assembled));
1976     CeedCall(CeedElemRestrictionReferenceCopy(assembled_rstr, rstr));
1977     CeedCall(CeedVectorDestroy(&assembled_vec));
1978     CeedCall(CeedElemRestrictionDestroy(&assembled_rstr));
1979   } else {
1980     // Operator fallback
1981     CeedOperator op_fallback;
1982 
1983     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
1984     if (op_fallback) CeedCall(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op_fallback, assembled, rstr, request));
1985     else return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedOperatorLinearAssembleQFunctionUpdate");
1986   }
1987   return CEED_ERROR_SUCCESS;
1988 }
1989 
1990 /**
1991   @brief Assemble the diagonal of a square linear `CeedOperator`
1992 
1993   This overwrites a `CeedVector` with the diagonal of a linear `CeedOperator`.
1994 
1995   Note: Currently only non-composite `CeedOperator` with a single field and composite `CeedOperator` with single field sub-operators are supported.
1996 
1997   Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
1998 
1999   @param[in]  op        `CeedOperator` to assemble `CeedQFunction`
2000   @param[out] assembled `CeedVector` to store assembled `CeedOperator` diagonal
2001   @param[in]  request   Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
2002 
2003   @return An error code: 0 - success, otherwise - failure
2004 
2005   @ref User
2006 **/
2007 int CeedOperatorLinearAssembleDiagonal(CeedOperator op, CeedVector assembled, CeedRequest *request) {
2008   bool     is_composite;
2009   CeedSize input_size = 0, output_size = 0;
2010 
2011   CeedCall(CeedOperatorCheckReady(op));
2012   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2013 
2014   CeedCall(CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
2015   CeedCheck(input_size == output_size, CeedOperatorReturnCeed(op), CEED_ERROR_DIMENSION, "Operator must be square");
2016 
2017   // Early exit for empty operator
2018   if (!is_composite) {
2019     CeedInt num_elem = 0;
2020 
2021     CeedCall(CeedOperatorGetNumElements(op, &num_elem));
2022     if (num_elem == 0) return CEED_ERROR_SUCCESS;
2023   }
2024 
2025   if (op->LinearAssembleDiagonal) {
2026     // Backend version
2027     CeedCall(op->LinearAssembleDiagonal(op, assembled, request));
2028     return CEED_ERROR_SUCCESS;
2029   } else if (op->LinearAssembleAddDiagonal) {
2030     // Backend version with zeroing first
2031     CeedCall(CeedVectorSetValue(assembled, 0.0));
2032     CeedCall(op->LinearAssembleAddDiagonal(op, assembled, request));
2033     return CEED_ERROR_SUCCESS;
2034   } else {
2035     // Operator fallback
2036     CeedOperator op_fallback;
2037 
2038     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2039     if (op_fallback) {
2040       CeedCall(CeedOperatorLinearAssembleDiagonal(op_fallback, assembled, request));
2041       return CEED_ERROR_SUCCESS;
2042     }
2043   }
2044   // Default interface implementation
2045   CeedCall(CeedVectorSetValue(assembled, 0.0));
2046   CeedCall(CeedOperatorLinearAssembleAddDiagonal(op, assembled, request));
2047   return CEED_ERROR_SUCCESS;
2048 }
2049 
2050 /**
2051   @brief Assemble the diagonal of a square linear `CeedOperator`.
2052 
2053   This sums into a `CeedVector` the diagonal of a linear `CeedOperator`.
2054 
2055   Note: Currently only non-composite `CeedOperator` with a single field and composite `CeedOperator` with single field sub-operators are supported.
2056 
2057   Note: Calling this function asserts that setup is complete and sets the CeedOperator as immutable.
2058 
2059   @param[in]  op        `CeedOperator` to assemble `CeedQFunction`
2060   @param[out] assembled `CeedVector` to store assembled `CeedOperator` diagonal
2061   @param[in]  request   Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
2062 
2063   @return An error code: 0 - success, otherwise - failure
2064 
2065   @ref User
2066 **/
2067 int CeedOperatorLinearAssembleAddDiagonal(CeedOperator op, CeedVector assembled, CeedRequest *request) {
2068   bool     is_composite;
2069   CeedSize input_size = 0, output_size = 0;
2070 
2071   CeedCall(CeedOperatorCheckReady(op));
2072   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2073 
2074   CeedCall(CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
2075   CeedCheck(input_size == output_size, CeedOperatorReturnCeed(op), CEED_ERROR_DIMENSION, "Operator must be square");
2076 
2077   // Early exit for empty operator
2078   if (!is_composite) {
2079     CeedInt num_elem = 0;
2080 
2081     CeedCall(CeedOperatorGetNumElements(op, &num_elem));
2082     if (num_elem == 0) return CEED_ERROR_SUCCESS;
2083   }
2084 
2085   if (op->LinearAssembleAddDiagonal) {
2086     // Backend version
2087     CeedCall(op->LinearAssembleAddDiagonal(op, assembled, request));
2088     return CEED_ERROR_SUCCESS;
2089   } else {
2090     // Operator fallback
2091     CeedOperator op_fallback;
2092 
2093     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2094     if (op_fallback) {
2095       CeedCall(CeedOperatorLinearAssembleAddDiagonal(op_fallback, assembled, request));
2096       return CEED_ERROR_SUCCESS;
2097     }
2098   }
2099   // Default interface implementation
2100   if (is_composite) {
2101     CeedCall(CeedCompositeOperatorLinearAssembleAddDiagonal(op, request, false, assembled));
2102   } else {
2103     CeedCall(CeedSingleOperatorLinearAssembleAddDiagonal(op, request, false, assembled));
2104   }
2105   return CEED_ERROR_SUCCESS;
2106 }
2107 
2108 /**
2109    @brief Fully assemble the point-block diagonal pattern of a linear `CeedOperator`.
2110 
2111    Expected to be used in conjunction with @ref CeedOperatorLinearAssemblePointBlockDiagonal().
2112 
2113    The assembly routines use coordinate format, with `num_entries` tuples of the form `(i, j, value)` which indicate that value should be added to the matrix in entry `(i, j)`.
2114    Note that the `(i, j)` pairs are unique.
2115    This function returns the number of entries and their `(i, j)` locations, while @ref CeedOperatorLinearAssemblePointBlockDiagonal() provides the values in the same ordering.
2116 
2117    This will generally be slow unless your operator is low-order.
2118 
2119    Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2120 
2121    @param[in]  op          `CeedOperator` to assemble
2122    @param[out] num_entries Number of entries in coordinate nonzero pattern
2123    @param[out] rows        Row number for each entry
2124    @param[out] cols        Column number for each entry
2125 
2126    @ref User
2127 **/
2128 int CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(CeedOperator op, CeedSize *num_entries, CeedInt **rows, CeedInt **cols) {
2129   bool          is_composite;
2130   CeedInt       num_active_components, num_sub_operators;
2131   CeedOperator *sub_operators;
2132 
2133   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2134 
2135   CeedSize input_size = 0, output_size = 0;
2136   CeedCall(CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
2137   CeedCheck(input_size == output_size, CeedOperatorReturnCeed(op), CEED_ERROR_DIMENSION, "Operator must be square");
2138 
2139   if (is_composite) {
2140     CeedCall(CeedCompositeOperatorGetNumSub(op, &num_sub_operators));
2141     CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
2142   } else {
2143     sub_operators     = &op;
2144     num_sub_operators = 1;
2145   }
2146 
2147   // Verify operator can be assembled correctly
2148   {
2149     CeedOperatorAssemblyData data;
2150     CeedInt                  num_active_elem_rstrs, comp_stride;
2151     CeedElemRestriction     *active_elem_rstrs;
2152 
2153     // Get initial values to check against
2154     CeedCall(CeedOperatorGetOperatorAssemblyData(sub_operators[0], &data));
2155     CeedCall(CeedOperatorAssemblyDataGetElemRestrictions(data, &num_active_elem_rstrs, &active_elem_rstrs, NULL, NULL));
2156     CeedCall(CeedElemRestrictionGetCompStride(active_elem_rstrs[0], &comp_stride));
2157     CeedCall(CeedElemRestrictionGetNumComponents(active_elem_rstrs[0], &num_active_components));
2158 
2159     // Verify that all active element restrictions have same component stride and number of components
2160     for (CeedInt k = 0; k < num_sub_operators; k++) {
2161       CeedCall(CeedOperatorGetOperatorAssemblyData(sub_operators[k], &data));
2162       CeedCall(CeedOperatorAssemblyDataGetElemRestrictions(data, &num_active_elem_rstrs, &active_elem_rstrs, NULL, NULL));
2163       for (CeedInt i = 0; i < num_active_elem_rstrs; i++) {
2164         CeedInt comp_stride_sub, num_active_components_sub;
2165 
2166         CeedCall(CeedElemRestrictionGetCompStride(active_elem_rstrs[i], &comp_stride_sub));
2167         CeedCheck(comp_stride == comp_stride_sub, CeedOperatorReturnCeed(op), CEED_ERROR_DIMENSION,
2168                   "Active element restrictions must have the same component stride: %d vs %d", comp_stride, comp_stride_sub);
2169         CeedCall(CeedElemRestrictionGetNumComponents(active_elem_rstrs[i], &num_active_components_sub));
2170         CeedCheck(num_active_components == num_active_components_sub, CeedOperatorReturnCeed(op), CEED_ERROR_INCOMPATIBLE,
2171                   "All suboperators must have the same number of output components."
2172                   " Previous: %" CeedInt_FMT " Current: %" CeedInt_FMT,
2173                   num_active_components, num_active_components_sub);
2174       }
2175     }
2176   }
2177   *num_entries = input_size * num_active_components;
2178   CeedCall(CeedCalloc(*num_entries, rows));
2179   CeedCall(CeedCalloc(*num_entries, cols));
2180 
2181   for (CeedInt o = 0; o < num_sub_operators; o++) {
2182     CeedElemRestriction active_elem_rstr, point_block_active_elem_rstr;
2183     CeedInt             comp_stride, num_elem, elem_size;
2184     const CeedInt      *offsets, *point_block_offsets;
2185 
2186     CeedCall(CeedOperatorGetActiveElemRestriction(sub_operators[o], &active_elem_rstr));
2187     CeedCall(CeedElemRestrictionGetCompStride(active_elem_rstr, &comp_stride));
2188     CeedCall(CeedElemRestrictionGetNumElements(active_elem_rstr, &num_elem));
2189     CeedCall(CeedElemRestrictionGetElementSize(active_elem_rstr, &elem_size));
2190     CeedCall(CeedElemRestrictionGetOffsets(active_elem_rstr, CEED_MEM_HOST, &offsets));
2191 
2192     CeedCall(CeedOperatorCreateActivePointBlockRestriction(active_elem_rstr, &point_block_active_elem_rstr));
2193     CeedCall(CeedElemRestrictionGetOffsets(point_block_active_elem_rstr, CEED_MEM_HOST, &point_block_offsets));
2194 
2195     for (CeedSize i = 0; i < num_elem * elem_size; i++) {
2196       for (CeedInt c_out = 0; c_out < num_active_components; c_out++) {
2197         for (CeedInt c_in = 0; c_in < num_active_components; c_in++) {
2198           (*rows)[point_block_offsets[i] + c_out * num_active_components + c_in] = offsets[i] + c_out * comp_stride;
2199           (*cols)[point_block_offsets[i] + c_out * num_active_components + c_in] = offsets[i] + c_in * comp_stride;
2200         }
2201       }
2202     }
2203 
2204     CeedCall(CeedElemRestrictionRestoreOffsets(active_elem_rstr, &offsets));
2205     CeedCall(CeedElemRestrictionRestoreOffsets(point_block_active_elem_rstr, &point_block_offsets));
2206     CeedCall(CeedElemRestrictionDestroy(&active_elem_rstr));
2207     CeedCall(CeedElemRestrictionDestroy(&point_block_active_elem_rstr));
2208   }
2209   return CEED_ERROR_SUCCESS;
2210 }
2211 
2212 /**
2213   @brief Assemble the point block diagonal of a square linear `CeedOperator`.
2214 
2215   This overwrites a `CeedVector` with the point block diagonal of a linear `CeedOperator`.
2216 
2217   Note: Currently only non-composite `CeedOperator` with a single field and composite `CeedOperator` with single field sub-operators are supported.
2218 
2219   Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2220 
2221   @param[in]  op        `CeedOperator` to assemble `CeedQFunction`
2222   @param[out] assembled `CeedVector` to store assembled `CeedOperator` point block diagonal, provided in row-major form with an `num_comp * num_comp` block at each node.
2223                           The dimensions of this vector are derived from the active vector for the `CeedOperator`.
2224                           The array has shape `[nodes, component out, component in]`.
2225   @param[in]  request   Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
2226 
2227   @return An error code: 0 - success, otherwise - failure
2228 
2229   @ref User
2230 **/
2231 int CeedOperatorLinearAssemblePointBlockDiagonal(CeedOperator op, CeedVector assembled, CeedRequest *request) {
2232   bool     is_composite;
2233   CeedSize input_size = 0, output_size = 0;
2234 
2235   CeedCall(CeedOperatorCheckReady(op));
2236   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2237 
2238   CeedCall(CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
2239   CeedCheck(input_size == output_size, CeedOperatorReturnCeed(op), CEED_ERROR_DIMENSION, "Operator must be square");
2240 
2241   // Early exit for empty operator
2242   if (!is_composite) {
2243     CeedInt num_elem = 0;
2244 
2245     CeedCall(CeedOperatorGetNumElements(op, &num_elem));
2246     if (num_elem == 0) return CEED_ERROR_SUCCESS;
2247   }
2248 
2249   if (op->LinearAssemblePointBlockDiagonal) {
2250     // Backend version
2251     CeedCall(op->LinearAssemblePointBlockDiagonal(op, assembled, request));
2252     return CEED_ERROR_SUCCESS;
2253   } else if (op->LinearAssembleAddPointBlockDiagonal) {
2254     // Backend version with zeroing first
2255     CeedCall(CeedVectorSetValue(assembled, 0.0));
2256     CeedCall(CeedOperatorLinearAssembleAddPointBlockDiagonal(op, assembled, request));
2257     return CEED_ERROR_SUCCESS;
2258   } else {
2259     // Operator fallback
2260     CeedOperator op_fallback;
2261 
2262     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2263     if (op_fallback) {
2264       CeedCall(CeedOperatorLinearAssemblePointBlockDiagonal(op_fallback, assembled, request));
2265       return CEED_ERROR_SUCCESS;
2266     }
2267   }
2268   // Default interface implementation
2269   CeedCall(CeedVectorSetValue(assembled, 0.0));
2270   CeedCall(CeedOperatorLinearAssembleAddPointBlockDiagonal(op, assembled, request));
2271   return CEED_ERROR_SUCCESS;
2272 }
2273 
2274 /**
2275   @brief Assemble the point block diagonal of a square linear `CeedOperator`.
2276 
2277   This sums into a `CeedVector` with the point block diagonal of a linear `CeedOperator`.
2278 
2279   Note: Currently only non-composite `CeedOperator` with a single field and composite `CeedOperator` with single field sub-operators are supported.
2280 
2281   Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2282 
2283   @param[in]  op        `CeedOperator` to assemble `CeedQFunction`
2284   @param[out] assembled `CeedVector` to store assembled CeedOperator point block diagonal, provided in row-major form with an `num_comp * num_comp` block at each node.
2285                           The dimensions of this vector are derived from the active vector for the `CeedOperator`.
2286                           The array has shape `[nodes, component out, component in]`.
2287   @param[in]  request   Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
2288 
2289   @return An error code: 0 - success, otherwise - failure
2290 
2291   @ref User
2292 **/
2293 int CeedOperatorLinearAssembleAddPointBlockDiagonal(CeedOperator op, CeedVector assembled, CeedRequest *request) {
2294   bool     is_composite;
2295   CeedSize input_size = 0, output_size = 0;
2296 
2297   CeedCall(CeedOperatorCheckReady(op));
2298   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2299 
2300   CeedCall(CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
2301   CeedCheck(input_size == output_size, CeedOperatorReturnCeed(op), CEED_ERROR_DIMENSION, "Operator must be square");
2302 
2303   // Early exit for empty operator
2304   if (!is_composite) {
2305     CeedInt num_elem = 0;
2306 
2307     CeedCall(CeedOperatorGetNumElements(op, &num_elem));
2308     if (num_elem == 0) return CEED_ERROR_SUCCESS;
2309   }
2310 
2311   if (op->LinearAssembleAddPointBlockDiagonal) {
2312     // Backend version
2313     CeedCall(op->LinearAssembleAddPointBlockDiagonal(op, assembled, request));
2314     return CEED_ERROR_SUCCESS;
2315   } else {
2316     // Operator fallback
2317     CeedOperator op_fallback;
2318 
2319     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2320     if (op_fallback) {
2321       CeedCall(CeedOperatorLinearAssembleAddPointBlockDiagonal(op_fallback, assembled, request));
2322       return CEED_ERROR_SUCCESS;
2323     }
2324   }
2325   // Default interface implementation
2326   if (is_composite) {
2327     CeedCall(CeedCompositeOperatorLinearAssembleAddDiagonal(op, request, true, assembled));
2328   } else {
2329     CeedCall(CeedSingleOperatorLinearAssembleAddDiagonal(op, request, true, assembled));
2330   }
2331   return CEED_ERROR_SUCCESS;
2332 }
2333 
2334 /**
2335    @brief Fully assemble the nonzero pattern of a linear `CeedOperator`.
2336 
2337    Expected to be used in conjunction with @ref CeedOperatorLinearAssemble().
2338 
2339    The assembly routines use coordinate format, with `num_entries` tuples of the form `(i, j, value)` which indicate that value should be added to the matrix in entry `(i, j)`.
2340    Note that the `(i, j)` pairs are not unique and may repeat.
2341    This function returns the number of entries and their `(i, j)` locations, while @ref CeedOperatorLinearAssemble() provides the values in the same ordering.
2342 
2343    This will generally be slow unless your operator is low-order.
2344 
2345    Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2346 
2347    @param[in]  op          `CeedOperator` to assemble
2348    @param[out] num_entries Number of entries in coordinate nonzero pattern
2349    @param[out] rows        Row number for each entry
2350    @param[out] cols        Column number for each entry
2351 
2352    @ref User
2353 **/
2354 int CeedOperatorLinearAssembleSymbolic(CeedOperator op, CeedSize *num_entries, CeedInt **rows, CeedInt **cols) {
2355   bool          is_composite;
2356   CeedInt       num_suboperators, offset = 0;
2357   CeedSize      single_entries;
2358   CeedOperator *sub_operators;
2359 
2360   CeedCall(CeedOperatorCheckReady(op));
2361   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2362 
2363   if (op->LinearAssembleSymbolic) {
2364     // Backend version
2365     CeedCall(op->LinearAssembleSymbolic(op, num_entries, rows, cols));
2366     return CEED_ERROR_SUCCESS;
2367   } else {
2368     // Operator fallback
2369     CeedOperator op_fallback;
2370 
2371     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2372     if (op_fallback) {
2373       CeedCall(CeedOperatorLinearAssembleSymbolic(op_fallback, num_entries, rows, cols));
2374       return CEED_ERROR_SUCCESS;
2375     }
2376   }
2377 
2378   // Default interface implementation
2379 
2380   // Count entries and allocate rows, cols arrays
2381   *num_entries = 0;
2382   if (is_composite) {
2383     CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators));
2384     CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
2385     for (CeedInt k = 0; k < num_suboperators; ++k) {
2386       CeedCall(CeedSingleOperatorAssemblyCountEntries(sub_operators[k], &single_entries));
2387       *num_entries += single_entries;
2388     }
2389   } else {
2390     CeedCall(CeedSingleOperatorAssemblyCountEntries(op, &single_entries));
2391     *num_entries += single_entries;
2392   }
2393   CeedCall(CeedCalloc(*num_entries, rows));
2394   CeedCall(CeedCalloc(*num_entries, cols));
2395 
2396   // Assemble nonzero locations
2397   if (is_composite) {
2398     CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators));
2399     CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
2400     for (CeedInt k = 0; k < num_suboperators; ++k) {
2401       CeedCall(CeedSingleOperatorAssembleSymbolic(sub_operators[k], offset, *rows, *cols));
2402       CeedCall(CeedSingleOperatorAssemblyCountEntries(sub_operators[k], &single_entries));
2403       offset += single_entries;
2404     }
2405   } else {
2406     CeedCall(CeedSingleOperatorAssembleSymbolic(op, offset, *rows, *cols));
2407   }
2408   return CEED_ERROR_SUCCESS;
2409 }
2410 
2411 /**
2412    @brief Fully assemble the nonzero entries of a linear operator.
2413 
2414    Expected to be used in conjunction with @ref CeedOperatorLinearAssembleSymbolic().
2415 
2416    The assembly routines use coordinate format, with `num_entries` tuples of the form `(i, j, value)` which indicate that value should be added to the matrix in entry `(i, j)`.
2417    Note that the `(i, j)` pairs are not unique and may repeat.
2418    This function returns the values of the nonzero entries to be added, their `(i, j)` locations are provided by @ref CeedOperatorLinearAssembleSymbolic().
2419 
2420    This will generally be slow unless your operator is low-order.
2421 
2422    Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2423 
2424    @param[in]  op     `CeedOperator` to assemble
2425    @param[out] values Values to assemble into matrix
2426 
2427    @ref User
2428 **/
2429 int CeedOperatorLinearAssemble(CeedOperator op, CeedVector values) {
2430   bool          is_composite;
2431   CeedInt       num_suboperators, offset = 0;
2432   CeedSize      single_entries = 0;
2433   CeedOperator *sub_operators;
2434 
2435   CeedCall(CeedOperatorCheckReady(op));
2436   CeedCall(CeedOperatorIsComposite(op, &is_composite));
2437 
2438   // Early exit for empty operator
2439   if (!is_composite) {
2440     CeedInt num_elem = 0;
2441 
2442     CeedCall(CeedOperatorGetNumElements(op, &num_elem));
2443     if (num_elem == 0) return CEED_ERROR_SUCCESS;
2444   }
2445 
2446   if (op->LinearAssemble) {
2447     // Backend version
2448     CeedCall(op->LinearAssemble(op, values));
2449     return CEED_ERROR_SUCCESS;
2450   } else {
2451     // Operator fallback
2452     CeedOperator op_fallback;
2453 
2454     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2455     if (op_fallback) {
2456       CeedCall(CeedOperatorLinearAssemble(op_fallback, values));
2457       return CEED_ERROR_SUCCESS;
2458     }
2459   }
2460 
2461   // Default interface implementation
2462   CeedCall(CeedVectorSetValue(values, 0.0));
2463   if (is_composite) {
2464     CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators));
2465     CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
2466     for (CeedInt k = 0; k < num_suboperators; k++) {
2467       CeedCall(CeedSingleOperatorAssemble(sub_operators[k], offset, values));
2468       CeedCall(CeedSingleOperatorAssemblyCountEntries(sub_operators[k], &single_entries));
2469       offset += single_entries;
2470     }
2471   } else {
2472     CeedCall(CeedSingleOperatorAssemble(op, offset, values));
2473   }
2474   return CEED_ERROR_SUCCESS;
2475 }
2476 
2477 /**
2478   @brief Get the multiplicity of nodes across sub-operators in a composite `CeedOperator`.
2479 
2480   Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2481 
2482   @param[in]  op               Composite `CeedOperator`
2483   @param[in]  num_skip_indices Number of sub-operators to skip
2484   @param[in]  skip_indices     Array of indices of sub-operators to skip
2485   @param[out] mult             Vector to store multiplicity (of size `l_size` )
2486 
2487   @return An error code: 0 - success, otherwise - failure
2488 
2489   @ref User
2490 **/
2491 int CeedCompositeOperatorGetMultiplicity(CeedOperator op, CeedInt num_skip_indices, CeedInt *skip_indices, CeedVector mult) {
2492   Ceed                ceed;
2493   CeedInt             num_suboperators;
2494   CeedSize            l_vec_len;
2495   CeedScalar         *mult_array;
2496   CeedVector          ones_l_vec;
2497   CeedElemRestriction elem_rstr, mult_elem_rstr;
2498   CeedOperator       *sub_operators;
2499 
2500   CeedCall(CeedOperatorCheckReady(op));
2501 
2502   // Zero mult vector
2503   CeedCall(CeedVectorSetValue(mult, 0.0));
2504 
2505   // Get suboperators
2506   CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators));
2507   if (num_suboperators == 0) return CEED_ERROR_SUCCESS;
2508   CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators));
2509 
2510   // Work vector
2511   CeedCall(CeedVectorGetLength(mult, &l_vec_len));
2512   CeedCall(CeedOperatorGetCeed(op, &ceed));
2513   CeedCall(CeedVectorCreate(ceed, l_vec_len, &ones_l_vec));
2514   CeedCall(CeedDestroy(&ceed));
2515   CeedCall(CeedVectorSetValue(ones_l_vec, 1.0));
2516   CeedCall(CeedVectorGetArray(mult, CEED_MEM_HOST, &mult_array));
2517 
2518   // Compute multiplicity across suboperators
2519   for (CeedInt i = 0; i < num_suboperators; i++) {
2520     const CeedScalar *sub_mult_array;
2521     CeedVector        sub_mult_l_vec, ones_e_vec;
2522 
2523     // -- Check for suboperator to skip
2524     for (CeedInt j = 0; j < num_skip_indices; j++) {
2525       if (skip_indices[j] == i) continue;
2526     }
2527 
2528     // -- Sub operator multiplicity
2529     CeedCall(CeedOperatorGetActiveElemRestriction(sub_operators[i], &elem_rstr));
2530     CeedCall(CeedElemRestrictionCreateUnorientedCopy(elem_rstr, &mult_elem_rstr));
2531     CeedCall(CeedElemRestrictionDestroy(&elem_rstr));
2532     CeedCall(CeedElemRestrictionCreateVector(mult_elem_rstr, &sub_mult_l_vec, &ones_e_vec));
2533     CeedCall(CeedVectorSetValue(sub_mult_l_vec, 0.0));
2534     CeedCall(CeedElemRestrictionApply(mult_elem_rstr, CEED_NOTRANSPOSE, ones_l_vec, ones_e_vec, CEED_REQUEST_IMMEDIATE));
2535     CeedCall(CeedElemRestrictionApply(mult_elem_rstr, CEED_TRANSPOSE, ones_e_vec, sub_mult_l_vec, CEED_REQUEST_IMMEDIATE));
2536     CeedCall(CeedVectorGetArrayRead(sub_mult_l_vec, CEED_MEM_HOST, &sub_mult_array));
2537     // ---- Flag every node present in the current suboperator
2538     for (CeedSize j = 0; j < l_vec_len; j++) {
2539       if (sub_mult_array[j] > 0.0) mult_array[j] += 1.0;
2540     }
2541     CeedCall(CeedVectorRestoreArrayRead(sub_mult_l_vec, &sub_mult_array));
2542     CeedCall(CeedVectorDestroy(&sub_mult_l_vec));
2543     CeedCall(CeedVectorDestroy(&ones_e_vec));
2544     CeedCall(CeedElemRestrictionDestroy(&mult_elem_rstr));
2545   }
2546   CeedCall(CeedVectorRestoreArray(mult, &mult_array));
2547   CeedCall(CeedVectorDestroy(&ones_l_vec));
2548   return CEED_ERROR_SUCCESS;
2549 }
2550 
2551 /**
2552   @brief Create a multigrid coarse `CeedOperator` and level transfer `CeedOperator` for a `CeedOperator`, creating the prolongation basis from the fine and coarse grid interpolation.
2553 
2554   Note: Calling this function asserts that setup is complete and sets all four `CeedOperator` as immutable.
2555 
2556   @param[in]  op_fine      Fine grid `CeedOperator`
2557   @param[in]  p_mult_fine  L-vector multiplicity in parallel gather/scatter, or `NULL` if not creating prolongation/restriction `CeedOperator`
2558   @param[in]  rstr_coarse  Coarse grid `CeedElemRestriction`
2559   @param[in]  basis_coarse Coarse grid active vector `CeedBasis`
2560   @param[out] op_coarse    Coarse grid `CeedOperator`
2561   @param[out] op_prolong   Coarse to fine `CeedOperator`, or `NULL`
2562   @param[out] op_restrict  Fine to coarse `CeedOperator`, or `NULL`
2563 
2564   @return An error code: 0 - success, otherwise - failure
2565 
2566   @ref User
2567 **/
2568 int CeedOperatorMultigridLevelCreate(CeedOperator op_fine, CeedVector p_mult_fine, CeedElemRestriction rstr_coarse, CeedBasis basis_coarse,
2569                                      CeedOperator *op_coarse, CeedOperator *op_prolong, CeedOperator *op_restrict) {
2570   CeedBasis basis_c_to_f = NULL;
2571 
2572   CeedCall(CeedOperatorCheckReady(op_fine));
2573 
2574   // Build prolongation matrix, if required
2575   if (op_prolong || op_restrict) {
2576     CeedBasis basis_fine;
2577 
2578     CeedCall(CeedOperatorGetActiveBasis(op_fine, &basis_fine));
2579     CeedCall(CeedBasisCreateProjection(basis_coarse, basis_fine, &basis_c_to_f));
2580     CeedCall(CeedBasisDestroy(&basis_fine));
2581   }
2582 
2583   // Core code
2584   CeedCall(CeedSingleOperatorMultigridLevel(op_fine, p_mult_fine, rstr_coarse, basis_coarse, basis_c_to_f, op_coarse, op_prolong, op_restrict));
2585   return CEED_ERROR_SUCCESS;
2586 }
2587 
2588 /**
2589   @brief Create a multigrid coarse `CeedOperator` and level transfer `CeedOperator` for a `CeedOperator` with a tensor basis for the active basis.
2590 
2591   Note: Calling this function asserts that setup is complete and sets all four `CeedOperator` as immutable.
2592 
2593   @param[in]  op_fine       Fine grid `CeedOperator`
2594   @param[in]  p_mult_fine   L-vector multiplicity in parallel gather/scatter, or `NULL` if not creating prolongation/restriction `CeedOperator`
2595   @param[in]  rstr_coarse   Coarse grid `CeedElemRestriction`
2596   @param[in]  basis_coarse  Coarse grid active vector `CeedBasis`
2597   @param[in]  interp_c_to_f Matrix for coarse to fine interpolation, or `NULL` if not creating prolongation/restriction `CeedOperator`
2598   @param[out] op_coarse     Coarse grid `CeedOperator`
2599   @param[out] op_prolong    Coarse to fine `CeedOperator`, or `NULL`
2600   @param[out] op_restrict   Fine to coarse `CeedOperator`, or `NULL`
2601 
2602   @return An error code: 0 - success, otherwise - failure
2603 
2604   @ref User
2605 **/
2606 int CeedOperatorMultigridLevelCreateTensorH1(CeedOperator op_fine, CeedVector p_mult_fine, CeedElemRestriction rstr_coarse, CeedBasis basis_coarse,
2607                                              const CeedScalar *interp_c_to_f, CeedOperator *op_coarse, CeedOperator *op_prolong,
2608                                              CeedOperator *op_restrict) {
2609   Ceed      ceed;
2610   CeedInt   Q_f, Q_c;
2611   CeedBasis basis_fine, basis_c_to_f = NULL;
2612 
2613   CeedCall(CeedOperatorCheckReady(op_fine));
2614   CeedCall(CeedOperatorGetCeed(op_fine, &ceed));
2615 
2616   // Check for compatible quadrature spaces
2617   CeedCall(CeedOperatorGetActiveBasis(op_fine, &basis_fine));
2618   CeedCall(CeedBasisGetNumQuadraturePoints(basis_fine, &Q_f));
2619   CeedCall(CeedBasisGetNumQuadraturePoints(basis_coarse, &Q_c));
2620   CeedCheck(Q_f == Q_c, ceed, CEED_ERROR_DIMENSION,
2621             "Bases must have compatible quadrature spaces."
2622             " Fine grid: %" CeedInt_FMT " points, Coarse grid: %" CeedInt_FMT " points",
2623             Q_f, Q_c);
2624 
2625   // Create coarse to fine basis, if required
2626   if (op_prolong || op_restrict) {
2627     CeedInt     dim, num_comp, num_nodes_c, P_1d_f, P_1d_c;
2628     CeedScalar *q_ref, *q_weight, *grad;
2629 
2630     // Check if interpolation matrix is provided
2631     CeedCheck(interp_c_to_f, ceed, CEED_ERROR_INCOMPATIBLE,
2632               "Prolongation or restriction operator creation requires coarse-to-fine interpolation matrix");
2633     CeedCall(CeedBasisGetDimension(basis_fine, &dim));
2634     CeedCall(CeedBasisGetNumComponents(basis_fine, &num_comp));
2635     CeedCall(CeedBasisGetNumNodes1D(basis_fine, &P_1d_f));
2636     CeedCall(CeedBasisDestroy(&basis_fine));
2637     CeedCall(CeedElemRestrictionGetElementSize(rstr_coarse, &num_nodes_c));
2638     P_1d_c = dim == 1 ? num_nodes_c : dim == 2 ? sqrt(num_nodes_c) : cbrt(num_nodes_c);
2639     CeedCall(CeedCalloc(P_1d_f, &q_ref));
2640     CeedCall(CeedCalloc(P_1d_f, &q_weight));
2641     CeedCall(CeedCalloc(P_1d_f * P_1d_c * dim, &grad));
2642     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_c, P_1d_f, interp_c_to_f, grad, q_ref, q_weight, &basis_c_to_f));
2643     CeedCall(CeedFree(&q_ref));
2644     CeedCall(CeedFree(&q_weight));
2645     CeedCall(CeedFree(&grad));
2646   }
2647 
2648   // Core code
2649   CeedCall(CeedSingleOperatorMultigridLevel(op_fine, p_mult_fine, rstr_coarse, basis_coarse, basis_c_to_f, op_coarse, op_prolong, op_restrict));
2650   CeedCall(CeedDestroy(&ceed));
2651   return CEED_ERROR_SUCCESS;
2652 }
2653 
2654 /**
2655   @brief Create a multigrid coarse `CeedOperator` and level transfer `CeedOperator` for a `CeedOperator` with a non-tensor basis for the active vector
2656 
2657   Note: Calling this function asserts that setup is complete and sets all four `CeedOperator` as immutable.
2658 
2659   @param[in]  op_fine       Fine grid `CeedOperator`
2660   @param[in]  p_mult_fine   L-vector multiplicity in parallel gather/scatter, or `NULL` if not creating prolongation/restriction `CeedOperator`
2661   @param[in]  rstr_coarse   Coarse grid `CeedElemRestriction`
2662   @param[in]  basis_coarse  Coarse grid active vector `CeedBasis`
2663   @param[in]  interp_c_to_f Matrix for coarse to fine interpolation, or `NULL` if not creating prolongation/restriction `CeedOperator`
2664   @param[out] op_coarse     Coarse grid `CeedOperator`
2665   @param[out] op_prolong    Coarse to fine `CeedOperator`, or `NULL`
2666   @param[out] op_restrict   Fine to coarse `CeedOperator`, or `NULL`
2667 
2668   @return An error code: 0 - success, otherwise - failure
2669 
2670   @ref User
2671 **/
2672 int CeedOperatorMultigridLevelCreateH1(CeedOperator op_fine, CeedVector p_mult_fine, CeedElemRestriction rstr_coarse, CeedBasis basis_coarse,
2673                                        const CeedScalar *interp_c_to_f, CeedOperator *op_coarse, CeedOperator *op_prolong,
2674                                        CeedOperator *op_restrict) {
2675   Ceed      ceed;
2676   CeedInt   Q_f, Q_c;
2677   CeedBasis basis_fine, basis_c_to_f = NULL;
2678 
2679   CeedCall(CeedOperatorCheckReady(op_fine));
2680   CeedCall(CeedOperatorGetCeed(op_fine, &ceed));
2681 
2682   // Check for compatible quadrature spaces
2683   CeedCall(CeedOperatorGetActiveBasis(op_fine, &basis_fine));
2684   CeedCall(CeedBasisGetNumQuadraturePoints(basis_fine, &Q_f));
2685   CeedCall(CeedBasisGetNumQuadraturePoints(basis_coarse, &Q_c));
2686   CeedCheck(Q_f == Q_c, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
2687 
2688   // Coarse to fine basis
2689   if (op_prolong || op_restrict) {
2690     CeedInt          dim, num_comp, num_nodes_c, num_nodes_f;
2691     CeedScalar      *q_ref, *q_weight, *grad;
2692     CeedElemTopology topo;
2693 
2694     // Check if interpolation matrix is provided
2695     CeedCheck(interp_c_to_f, ceed, CEED_ERROR_INCOMPATIBLE,
2696               "Prolongation or restriction operator creation requires coarse-to-fine interpolation matrix");
2697     CeedCall(CeedBasisGetTopology(basis_fine, &topo));
2698     CeedCall(CeedBasisGetDimension(basis_fine, &dim));
2699     CeedCall(CeedBasisGetNumComponents(basis_fine, &num_comp));
2700     CeedCall(CeedBasisGetNumNodes(basis_fine, &num_nodes_f));
2701     CeedCall(CeedBasisDestroy(&basis_fine));
2702     CeedCall(CeedElemRestrictionGetElementSize(rstr_coarse, &num_nodes_c));
2703     CeedCall(CeedCalloc(num_nodes_f * dim, &q_ref));
2704     CeedCall(CeedCalloc(num_nodes_f, &q_weight));
2705     CeedCall(CeedCalloc(num_nodes_f * num_nodes_c * dim, &grad));
2706     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_c, num_nodes_f, interp_c_to_f, grad, q_ref, q_weight, &basis_c_to_f));
2707     CeedCall(CeedFree(&q_ref));
2708     CeedCall(CeedFree(&q_weight));
2709     CeedCall(CeedFree(&grad));
2710   }
2711 
2712   // Core code
2713   CeedCall(CeedSingleOperatorMultigridLevel(op_fine, p_mult_fine, rstr_coarse, basis_coarse, basis_c_to_f, op_coarse, op_prolong, op_restrict));
2714   CeedCall(CeedDestroy(&ceed));
2715   return CEED_ERROR_SUCCESS;
2716 }
2717 
2718 /**
2719   @brief Build a FDM based approximate inverse for each element for a `CeedOperator`.
2720 
2721   This returns a `CeedOperator` and `CeedVector` to apply a Fast Diagonalization Method based approximate inverse.
2722   This function obtains the simultaneous diagonalization for the 1D mass and Laplacian operators, \f$M = V^T V, K = V^T S V\f$.
2723   The assembled `CeedQFunction` is used to modify the eigenvalues from simultaneous diagonalization and obtain an approximate inverse of the form \f$V^T \hat S V\f$.
2724   The `CeedOperator` must be linear and non-composite.
2725   The associated `CeedQFunction` must therefore also be linear.
2726 
2727   Note: Calling this function asserts that setup is complete and sets the `CeedOperator` as immutable.
2728 
2729   @param[in]  op      `CeedOperator` to create element inverses
2730   @param[out] fdm_inv `CeedOperator` to apply the action of a FDM based inverse for each element
2731   @param[in]  request Address of @ref CeedRequest for non-blocking completion, else @ref CEED_REQUEST_IMMEDIATE
2732 
2733   @return An error code: 0 - success, otherwise - failure
2734 
2735   @ref User
2736 **/
2737 int CeedOperatorCreateFDMElementInverse(CeedOperator op, CeedOperator *fdm_inv, CeedRequest *request) {
2738   Ceed                 ceed, ceed_parent;
2739   bool                 interp = false, grad = false, is_tensor_basis = true;
2740   CeedInt              num_input_fields, P_1d, Q_1d, num_nodes, num_qpts, dim, num_comp = 1, num_elem = 1;
2741   CeedScalar          *mass, *laplace, *x, *fdm_interp, *lambda, *elem_avg;
2742   const CeedScalar    *interp_1d, *grad_1d, *q_weight_1d;
2743   CeedVector           q_data;
2744   CeedElemRestriction  rstr  = NULL, rstr_qd_i;
2745   CeedBasis            basis = NULL, fdm_basis;
2746   CeedQFunctionContext ctx_fdm;
2747   CeedQFunctionField  *qf_fields;
2748   CeedQFunction        qf, qf_fdm;
2749   CeedOperatorField   *op_fields;
2750 
2751   CeedCall(CeedOperatorCheckReady(op));
2752 
2753   if (op->CreateFDMElementInverse) {
2754     // Backend version
2755     CeedCall(op->CreateFDMElementInverse(op, fdm_inv, request));
2756     return CEED_ERROR_SUCCESS;
2757   } else {
2758     // Operator fallback
2759     CeedOperator op_fallback;
2760 
2761     CeedCall(CeedOperatorGetFallback(op, &op_fallback));
2762     if (op_fallback) {
2763       CeedCall(CeedOperatorCreateFDMElementInverse(op_fallback, fdm_inv, request));
2764       return CEED_ERROR_SUCCESS;
2765     }
2766   }
2767 
2768   // Default interface implementation
2769   CeedCall(CeedOperatorGetCeed(op, &ceed));
2770   CeedCall(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
2771   CeedCall(CeedOperatorGetQFunction(op, &qf));
2772 
2773   // Determine active input basis
2774   CeedCall(CeedOperatorGetFields(op, &num_input_fields, &op_fields, NULL, NULL));
2775   CeedCall(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
2776   for (CeedInt i = 0; i < num_input_fields; i++) {
2777     CeedVector vec;
2778 
2779     CeedCall(CeedOperatorFieldGetVector(op_fields[i], &vec));
2780     if (vec == CEED_VECTOR_ACTIVE) {
2781       CeedEvalMode eval_mode;
2782 
2783       CeedCall(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
2784       interp = interp || eval_mode == CEED_EVAL_INTERP;
2785       grad   = grad || eval_mode == CEED_EVAL_GRAD;
2786       if (!basis) CeedCall(CeedOperatorFieldGetBasis(op_fields[i], &basis));
2787       if (!rstr) CeedCall(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
2788     }
2789     CeedCall(CeedVectorDestroy(&vec));
2790   }
2791   CeedCheck(basis, ceed, CEED_ERROR_BACKEND, "No active field set");
2792   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
2793   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
2794   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
2795   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
2796   CeedCall(CeedBasisGetDimension(basis, &dim));
2797   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
2798   CeedCall(CeedElemRestrictionGetNumElements(rstr, &num_elem));
2799 
2800   // Build and diagonalize 1D Mass and Laplacian
2801   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
2802   CeedCheck(is_tensor_basis, ceed, CEED_ERROR_BACKEND, "FDMElementInverse only supported for tensor bases");
2803   CeedCall(CeedCalloc(P_1d * P_1d, &mass));
2804   CeedCall(CeedCalloc(P_1d * P_1d, &laplace));
2805   CeedCall(CeedCalloc(P_1d * P_1d, &x));
2806   CeedCall(CeedCalloc(P_1d * P_1d, &fdm_interp));
2807   CeedCall(CeedCalloc(P_1d, &lambda));
2808   // -- Build matrices
2809   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
2810   CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
2811   CeedCall(CeedBasisGetQWeights(basis, &q_weight_1d));
2812   CeedCall(CeedBuildMassLaplace(interp_1d, grad_1d, q_weight_1d, P_1d, Q_1d, dim, mass, laplace));
2813 
2814   // -- Diagonalize
2815   CeedCall(CeedSimultaneousDiagonalization(ceed, laplace, mass, x, lambda, P_1d));
2816   CeedCall(CeedFree(&mass));
2817   CeedCall(CeedFree(&laplace));
2818   for (CeedInt i = 0; i < P_1d; i++) {
2819     for (CeedInt j = 0; j < P_1d; j++) fdm_interp[i + j * P_1d] = x[j + i * P_1d];
2820   }
2821   CeedCall(CeedFree(&x));
2822 
2823   {
2824     CeedInt             layout[3], num_modes = (interp ? 1 : 0) + (grad ? dim : 0);
2825     CeedScalar          max_norm = 0;
2826     const CeedScalar   *assembled_array, *q_weight_array;
2827     CeedVector          assembled = NULL, q_weight;
2828     CeedElemRestriction rstr_qf   = NULL;
2829 
2830     // Assemble QFunction
2831     CeedCall(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled, &rstr_qf, request));
2832     CeedCall(CeedElemRestrictionGetELayout(rstr_qf, layout));
2833     CeedCall(CeedElemRestrictionDestroy(&rstr_qf));
2834     CeedCall(CeedVectorNorm(assembled, CEED_NORM_MAX, &max_norm));
2835 
2836     // Calculate element averages
2837     CeedCall(CeedVectorCreate(ceed_parent, num_qpts, &q_weight));
2838     CeedCall(CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_weight));
2839     CeedCall(CeedVectorGetArrayRead(assembled, CEED_MEM_HOST, &assembled_array));
2840     CeedCall(CeedVectorGetArrayRead(q_weight, CEED_MEM_HOST, &q_weight_array));
2841     CeedCall(CeedCalloc(num_elem, &elem_avg));
2842     const CeedScalar qf_value_bound = max_norm * 100 * CEED_EPSILON;
2843 
2844     for (CeedInt e = 0; e < num_elem; e++) {
2845       CeedInt count = 0;
2846 
2847       for (CeedInt q = 0; q < num_qpts; q++) {
2848         for (CeedInt i = 0; i < num_comp * num_comp * num_modes * num_modes; i++) {
2849           if (fabs(assembled_array[q * layout[0] + i * layout[1] + e * layout[2]]) > qf_value_bound) {
2850             elem_avg[e] += assembled_array[q * layout[0] + i * layout[1] + e * layout[2]] / q_weight_array[q];
2851             count++;
2852           }
2853         }
2854       }
2855       if (count) {
2856         elem_avg[e] /= count;
2857       } else {
2858         elem_avg[e] = 1.0;
2859       }
2860     }
2861     CeedCall(CeedVectorRestoreArrayRead(assembled, &assembled_array));
2862     CeedCall(CeedVectorDestroy(&assembled));
2863     CeedCall(CeedVectorRestoreArrayRead(q_weight, &q_weight_array));
2864     CeedCall(CeedVectorDestroy(&q_weight));
2865   }
2866 
2867   // Build FDM diagonal
2868   {
2869     CeedScalar *q_data_array, *fdm_diagonal;
2870 
2871     CeedCall(CeedCalloc(num_comp * num_nodes, &fdm_diagonal));
2872     const CeedScalar fdm_diagonal_bound = num_nodes * CEED_EPSILON;
2873     for (CeedInt c = 0; c < num_comp; c++) {
2874       for (CeedInt n = 0; n < num_nodes; n++) {
2875         if (interp) fdm_diagonal[c * num_nodes + n] = 1.0;
2876         if (grad) {
2877           for (CeedInt d = 0; d < dim; d++) {
2878             CeedInt i = (n / CeedIntPow(P_1d, d)) % P_1d;
2879             fdm_diagonal[c * num_nodes + n] += lambda[i];
2880           }
2881         }
2882         if (fabs(fdm_diagonal[c * num_nodes + n]) < fdm_diagonal_bound) fdm_diagonal[c * num_nodes + n] = fdm_diagonal_bound;
2883       }
2884     }
2885     CeedCall(CeedVectorCreate(ceed_parent, num_elem * num_comp * num_nodes, &q_data));
2886     CeedCall(CeedVectorSetValue(q_data, 0.0));
2887     CeedCall(CeedVectorGetArrayWrite(q_data, CEED_MEM_HOST, &q_data_array));
2888     for (CeedInt e = 0; e < num_elem; e++) {
2889       for (CeedInt c = 0; c < num_comp; c++) {
2890         for (CeedInt n = 0; n < num_nodes; n++) {
2891           q_data_array[(e * num_comp + c) * num_nodes + n] = 1. / (elem_avg[e] * fdm_diagonal[c * num_nodes + n]);
2892         }
2893       }
2894     }
2895     CeedCall(CeedFree(&elem_avg));
2896     CeedCall(CeedFree(&fdm_diagonal));
2897     CeedCall(CeedVectorRestoreArray(q_data, &q_data_array));
2898   }
2899 
2900   // Setup FDM operator
2901   // -- Basis
2902   {
2903     CeedScalar *grad_dummy, *q_ref_dummy, *q_weight_dummy;
2904 
2905     CeedCall(CeedCalloc(P_1d * P_1d, &grad_dummy));
2906     CeedCall(CeedCalloc(P_1d, &q_ref_dummy));
2907     CeedCall(CeedCalloc(P_1d, &q_weight_dummy));
2908     CeedCall(CeedBasisCreateTensorH1(ceed_parent, dim, num_comp, P_1d, P_1d, fdm_interp, grad_dummy, q_ref_dummy, q_weight_dummy, &fdm_basis));
2909     CeedCall(CeedFree(&fdm_interp));
2910     CeedCall(CeedFree(&grad_dummy));
2911     CeedCall(CeedFree(&q_ref_dummy));
2912     CeedCall(CeedFree(&q_weight_dummy));
2913     CeedCall(CeedFree(&lambda));
2914   }
2915 
2916   // -- Restriction
2917   {
2918     CeedInt strides[3] = {1, num_nodes, num_nodes * num_comp};
2919     CeedCall(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, num_nodes, num_comp,
2920                                               (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes, strides, &rstr_qd_i));
2921   }
2922 
2923   // -- QFunction
2924   CeedCall(CeedQFunctionCreateInteriorByName(ceed_parent, "Scale", &qf_fdm));
2925   CeedCall(CeedQFunctionAddInput(qf_fdm, "input", num_comp, CEED_EVAL_INTERP));
2926   CeedCall(CeedQFunctionAddInput(qf_fdm, "scale", num_comp, CEED_EVAL_NONE));
2927   CeedCall(CeedQFunctionAddOutput(qf_fdm, "output", num_comp, CEED_EVAL_INTERP));
2928   CeedCall(CeedQFunctionSetUserFlopsEstimate(qf_fdm, num_comp));
2929 
2930   // -- QFunction context
2931   {
2932     CeedInt *num_comp_data;
2933 
2934     CeedCall(CeedCalloc(1, &num_comp_data));
2935     num_comp_data[0] = num_comp;
2936     CeedCall(CeedQFunctionContextCreate(ceed, &ctx_fdm));
2937     CeedCall(CeedQFunctionContextSetData(ctx_fdm, CEED_MEM_HOST, CEED_OWN_POINTER, sizeof(*num_comp_data), num_comp_data));
2938   }
2939   CeedCall(CeedQFunctionSetContext(qf_fdm, ctx_fdm));
2940   CeedCall(CeedQFunctionContextDestroy(&ctx_fdm));
2941 
2942   // -- Operator
2943   CeedCall(CeedOperatorCreate(ceed_parent, qf_fdm, NULL, NULL, fdm_inv));
2944   CeedCall(CeedOperatorSetField(*fdm_inv, "input", rstr, fdm_basis, CEED_VECTOR_ACTIVE));
2945   CeedCall(CeedOperatorSetField(*fdm_inv, "scale", rstr_qd_i, CEED_BASIS_NONE, q_data));
2946   CeedCall(CeedOperatorSetField(*fdm_inv, "output", rstr, fdm_basis, CEED_VECTOR_ACTIVE));
2947 
2948   // Cleanup
2949   CeedCall(CeedDestroy(&ceed));
2950   CeedCall(CeedDestroy(&ceed_parent));
2951   CeedCall(CeedVectorDestroy(&q_data));
2952   CeedCall(CeedElemRestrictionDestroy(&rstr));
2953   CeedCall(CeedElemRestrictionDestroy(&rstr_qd_i));
2954   CeedCall(CeedBasisDestroy(&basis));
2955   CeedCall(CeedBasisDestroy(&fdm_basis));
2956   CeedCall(CeedQFunctionDestroy(&qf));
2957   CeedCall(CeedQFunctionDestroy(&qf_fdm));
2958   return CEED_ERROR_SUCCESS;
2959 }
2960 
2961 /// @}
2962