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