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