xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-qfunction.c (revision 9ba83ac0e4b1fca39d6fa6737a318a9f0cbc172d)
1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed-impl.h>
9 #include <ceed.h>
10 #include <ceed/backend.h>
11 #include <ceed/jit-tools.h>
12 #include <limits.h>
13 #include <stdbool.h>
14 #include <stdio.h>
15 #include <string.h>
16 
17 /// @file
18 /// Implementation of public CeedQFunction interfaces
19 
20 /// @cond DOXYGEN_SKIP
21 static struct CeedQFunction_private ceed_qfunction_none;
22 /// @endcond
23 
24 /// @addtogroup CeedQFunctionUser
25 /// @{
26 
27 // Indicate that no `CeedQFunction` is provided by the user
28 const CeedQFunction CEED_QFUNCTION_NONE = &ceed_qfunction_none;
29 
30 /// @}
31 
32 /// @cond DOXYGEN_SKIP
33 static struct {
34   char              name[CEED_MAX_RESOURCE_LEN];
35   char              source[CEED_MAX_RESOURCE_LEN];
36   CeedInt           vec_length;
37   CeedQFunctionUser f;
38   int (*init)(Ceed ceed, const char *name, CeedQFunction qf);
39 } gallery_qfunctions[1024];
40 static size_t num_qfunctions;
41 /// @endcond
42 
43 /// ----------------------------------------------------------------------------
44 /// CeedQFunction Library Internal Functions
45 /// ----------------------------------------------------------------------------
46 /// @addtogroup CeedQFunctionDeveloper
47 /// @{
48 
49 /**
50   @brief Register a gallery `CeedQFunction`
51 
52   @param[in] name       Name for this backend to respond to
53   @param[in] source     Absolute path to source of `CeedQFunction`, "\path\CEED_DIR\gallery\folder\file.h:function_name"
54   @param[in] vec_length Vector length.
55                           Caller must ensure that number of quadrature points is a multiple of `vec_length`.
56   @param[in] f          Function pointer to evaluate action at quadrature points.
57                           See `CeedQFunctionUser`.
58   @param[in] init       Initialization function called by @ref CeedQFunctionCreateInteriorByName() when the `CeedQFunction` is selected.
59 
60   @return An error code: 0 - success, otherwise - failure
61 
62   @ref Developer
63 **/
64 int CeedQFunctionRegister(const char *name, const char *source, CeedInt vec_length, CeedQFunctionUser f,
65                           int (*init)(Ceed, const char *, CeedQFunction)) {
66   const char *relative_file_path;
67   int         ierr = 0;
68 
69   CeedDebugEnv("Gallery Register: %s", name);
70   CeedCall(CeedGetJitRelativePath(source, &relative_file_path));
71   CeedPragmaCritical(CeedQFunctionRegister) {
72     if (num_qfunctions < sizeof(gallery_qfunctions) / sizeof(gallery_qfunctions[0])) {
73       strncpy(gallery_qfunctions[num_qfunctions].name, name, CEED_MAX_RESOURCE_LEN);
74       gallery_qfunctions[num_qfunctions].name[CEED_MAX_RESOURCE_LEN - 1] = 0;
75       strncpy(gallery_qfunctions[num_qfunctions].source, relative_file_path, CEED_MAX_RESOURCE_LEN);
76       gallery_qfunctions[num_qfunctions].source[CEED_MAX_RESOURCE_LEN - 1] = 0;
77       gallery_qfunctions[num_qfunctions].vec_length                        = vec_length;
78       gallery_qfunctions[num_qfunctions].f                                 = f;
79       gallery_qfunctions[num_qfunctions].init                              = init;
80       num_qfunctions++;
81     } else {
82       ierr = 1;
83     }
84   }
85   CeedCheck(ierr == 0, NULL, CEED_ERROR_MAJOR, "Too many gallery CeedQFunctions");
86   return CEED_ERROR_SUCCESS;
87 }
88 
89 /**
90   @brief Set a `CeedQFunction` field, used by @ref CeedQFunctionAddInput() and @ref CeedQFunctionAddOutput()
91 
92   @param[out] f           `CeedQFunctionField`
93   @param[in]  field_name  Name of `CeedQFunction` field
94   @param[in]  size        Size of `CeedQFunction` field, (`num_comp * 1`) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_WEIGHT, (`num_comp * 1`) for @ref CEED_EVAL_INTERP for an \f$H^1\f$ space or (`num_comp * dim`) for an \f$H(\mathrm{div})\f$ or \f$H(\mathrm{curl})\f$ space, (`num_comp * dim`) for @ref CEED_EVAL_GRAD, or (num_comp * 1) for @ref CEED_EVAL_DIV, and (`num_comp * curl_dim`) with `curl_dim = 1` if `dim < 3` and `curl_dim = dim` for @ref CEED_EVAL_CURL.
95   @param[in]  eval_mode   @ref CEED_EVAL_NONE to use values directly,
96                             @ref CEED_EVAL_WEIGHT to use quadrature weights,
97                             @ref CEED_EVAL_INTERP to use interpolated values,
98                             @ref CEED_EVAL_GRAD to use gradients,
99                             @ref CEED_EVAL_DIV to use divergence,
100                             @ref CEED_EVAL_CURL to use curl
101 
102   @return An error code: 0 - success, otherwise - failure
103 
104   @ref Developer
105 **/
106 static int CeedQFunctionFieldSet(CeedQFunctionField *f, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
107   CeedCall(CeedCalloc(1, f));
108   CeedCall(CeedStringAllocCopy(field_name, (char **)&(*f)->field_name));
109   (*f)->size      = size;
110   (*f)->eval_mode = eval_mode;
111   return CEED_ERROR_SUCCESS;
112 }
113 
114 /**
115   @brief View a field of a `CeedQFunction`
116 
117   @param[in] field        `CeedQFunction` field to view
118   @param[in] field_number Number of field being viewed
119   @param[in] in           true for input field, false for output
120   @param[in] tabs         Tabs to append before each new line
121   @param[in] stream       Stream to view to, e.g., `stdout`
122 
123   @return An error code: 0 - success, otherwise - failure
124 
125   @ref Utility
126 **/
127 static int CeedQFunctionFieldView(CeedQFunctionField field, CeedInt field_number, bool in, const char *tabs, FILE *stream) {
128   const char  *inout = in ? "Input" : "Output";
129   const char  *field_name;
130   CeedInt      size;
131   CeedEvalMode eval_mode;
132 
133   CeedCall(CeedQFunctionFieldGetData(field, &field_name, &size, &eval_mode));
134   fprintf(stream,
135           "%s    %s field %" CeedInt_FMT
136           ":\n%s"
137           "      Name: \"%s\"\n%s"
138           "      Size: %" CeedInt_FMT
139           "\n%s"
140           "      EvalMode: \"%s\"\n",
141           tabs, inout, field_number, tabs, field_name, tabs, size, tabs, CeedEvalModes[eval_mode]);
142   return CEED_ERROR_SUCCESS;
143 }
144 
145 /**
146   @brief Set flag to determine if Fortran interface is used
147 
148   @param[in,out] qf     CeedQFunction
149   @param[in]     status Boolean value to set as Fortran status
150 
151   @return An error code: 0 - success, otherwise - failure
152 
153   @ref Backend
154 **/
155 int CeedQFunctionSetFortranStatus(CeedQFunction qf, bool status) {
156   qf->is_fortran = status;
157   return CEED_ERROR_SUCCESS;
158 }
159 
160 /// @}
161 
162 /// ----------------------------------------------------------------------------
163 /// CeedQFunction Backend API
164 /// ----------------------------------------------------------------------------
165 /// @addtogroup CeedQFunctionBackend
166 /// @{
167 
168 /**
169   @brief Get the vector length of a `CeedQFunction`
170 
171   @param[in]  qf         `CeedQFunction`
172   @param[out] vec_length Variable to store vector length
173 
174   @return An error code: 0 - success, otherwise - failure
175 
176   @ref Backend
177 **/
178 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vec_length) {
179   *vec_length = qf->vec_length;
180   return CEED_ERROR_SUCCESS;
181 }
182 
183 /**
184   @brief Get the number of inputs and outputs to a `CeedQFunction`
185 
186   @param[in]  qf         `CeedQFunction`
187   @param[out] num_input  Variable to store number of input fields
188   @param[out] num_output Variable to store number of output fields
189 
190   @return An error code: 0 - success, otherwise - failure
191 
192   @ref Backend
193 **/
194 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *num_input, CeedInt *num_output) {
195   if (num_input) *num_input = qf->num_input_fields;
196   if (num_output) *num_output = qf->num_output_fields;
197   return CEED_ERROR_SUCCESS;
198 }
199 
200 /**
201   @brief Get the name of the `CeedQFunction`.
202     Use the `name` if created via @ref `CeedQFunctionCreateInteriorByName`, otherwise return the kernel name via @ref `CeedQFunctionGetKernelName`.
203 
204   @param[in]  qf          `CeedQFunction`
205   @param[out] kernel_name Variable to store `CeedQFunction` name
206 
207   @return An error code: 0 - success, otherwise - failure
208 
209   @ref Backend
210 **/
211 int CeedQFunctionGetName(CeedQFunction qf, const char **name) {
212   if (qf->is_gallery) {
213     *name = qf->gallery_name;
214   } else {
215     CeedCall(CeedQFunctionGetKernelName(qf, name));
216   }
217   return CEED_ERROR_SUCCESS;
218 }
219 
220 /**
221   @brief Get the name of the user function for a `CeedQFunction`
222 
223   @param[in]  qf          `CeedQFunction`
224   @param[out] kernel_name Variable to store string holding kernel name
225 
226   @return An error code: 0 - success, otherwise - failure
227 
228   @ref Backend
229 **/
230 int CeedQFunctionGetKernelName(CeedQFunction qf, const char **kernel_name) {
231   if (!qf->kernel_name) {
232     char *kernel_name_copy;
233 
234     if (qf->user_source) {
235       const char *kernel_name     = strrchr(qf->user_source, ':') + 1;
236       size_t      kernel_name_len = strlen(kernel_name);
237 
238       CeedCall(CeedCalloc(kernel_name_len + 1, &kernel_name_copy));
239       memcpy(kernel_name_copy, kernel_name, kernel_name_len);
240     } else {
241       CeedCall(CeedCalloc(1, &kernel_name_copy));
242     }
243     qf->kernel_name = kernel_name_copy;
244   }
245 
246   *kernel_name = qf->kernel_name;
247   return CEED_ERROR_SUCCESS;
248 }
249 
250 /**
251   @brief Get the source path string for a `CeedQFunction`
252 
253   @param[in]  qf          `CeedQFunction`
254   @param[out] source_path Variable to store source path string
255 
256   @return An error code: 0 - success, otherwise - failure
257 
258   @ref Backend
259 **/
260 int CeedQFunctionGetSourcePath(CeedQFunction qf, const char **source_path) {
261   if (!qf->source_path && qf->user_source) {
262     Ceed        ceed;
263     bool        is_absolute_path;
264     char       *source_path_copy;
265     const char *absolute_path;
266     const char *kernel_name     = strrchr(qf->user_source, ':') + 1;
267     size_t      kernel_name_len = strlen(kernel_name);
268 
269     CeedCall(CeedQFunctionGetCeed(qf, &ceed));
270     CeedCall(CeedCheckFilePath(ceed, qf->user_source, &is_absolute_path));
271     if (is_absolute_path) {
272       absolute_path = (char *)qf->user_source;
273     } else {
274       CeedCall(CeedGetJitAbsolutePath(ceed, qf->user_source, &absolute_path));
275     }
276     CeedCall(CeedDestroy(&ceed));
277 
278     size_t source_len = strlen(absolute_path) - kernel_name_len - 1;
279 
280     CeedCall(CeedCalloc(source_len + 1, &source_path_copy));
281     memcpy(source_path_copy, absolute_path, source_len);
282     qf->source_path = source_path_copy;
283 
284     if (!is_absolute_path) CeedCall(CeedFree(&absolute_path));
285   }
286 
287   *source_path = (char *)qf->source_path;
288   return CEED_ERROR_SUCCESS;
289 }
290 
291 /**
292   @brief Initialize and load `CeedQFunction` source file into string buffer, including full text of local files in place of `#include "local.h"`.
293 
294   The `buffer` is set to `NULL` if there is no `CeedQFunction` source file.
295 
296   Note: This function may as well return a mutable buffer, but all current uses
297   do not modify it. (This is just a downside of `const` semantics with output
298   arguments instead of returns.)
299 
300   Note: Caller is responsible for freeing the string buffer with @ref CeedFree().
301 
302   @param[in]  qf            `CeedQFunction`
303   @param[out] source_buffer String buffer for source file contents
304 
305   @return An error code: 0 - success, otherwise - failure
306 
307   @ref Backend
308 **/
309 int CeedQFunctionLoadSourceToBuffer(CeedQFunction qf, const char **source_buffer) {
310   const char *source_path;
311 
312   CeedCall(CeedQFunctionGetSourcePath(qf, &source_path));
313   *source_buffer = NULL;
314   if (source_path) {
315     Ceed  ceed;
316     char *buffer = NULL;
317 
318     CeedCall(CeedQFunctionGetCeed(qf, &ceed));
319     CeedCall(CeedLoadSourceToBuffer(ceed, source_path, &buffer));
320     CeedCall(CeedDestroy(&ceed));
321     *source_buffer = buffer;
322   }
323   return CEED_ERROR_SUCCESS;
324 }
325 
326 /**
327   @brief Get the User Function for a `CeedQFunction`
328 
329   @param[in]  qf `CeedQFunction`
330   @param[out] f  Variable to store user function
331 
332   @return An error code: 0 - success, otherwise - failure
333 
334   @ref Backend
335 **/
336 int CeedQFunctionGetUserFunction(CeedQFunction qf, CeedQFunctionUser *f) {
337   *f = qf->function;
338   return CEED_ERROR_SUCCESS;
339 }
340 
341 /**
342   @brief Get global context for a `CeedQFunction`.
343 
344   Note: For `CeedQFunction` from the Fortran interface, this function will return the Fortran context `CeedQFunctionContext`.
345 
346   @param[in]  qf  CeedQFunction
347   @param[out] ctx Variable to store CeedQFunctionContext
348 
349   @return An error code: 0 - success, otherwise - failure
350 
351   @ref Backend
352 **/
353 int CeedQFunctionGetContext(CeedQFunction qf, CeedQFunctionContext *ctx) {
354   *ctx = NULL;
355   if (qf->ctx) CeedCall(CeedQFunctionContextReferenceCopy(qf->ctx, ctx));
356   return CEED_ERROR_SUCCESS;
357 }
358 
359 /**
360   @brief Get context data of a `CeedQFunction`
361 
362   @param[in]  qf       `CeedQFunction`
363   @param[in]  mem_type Memory type on which to access the data.
364                          If the backend uses a different memory type, this will perform a copy.
365   @param[out] data     Data on memory type mem_type
366 
367   @return An error code: 0 - success, otherwise - failure
368 
369   @ref Backend
370 **/
371 int CeedQFunctionGetContextData(CeedQFunction qf, CeedMemType mem_type, void *data) {
372   bool                 is_writable;
373   CeedQFunctionContext ctx;
374 
375   CeedCall(CeedQFunctionGetContext(qf, &ctx));
376   if (ctx) {
377     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
378     if (is_writable) {
379       CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data));
380     } else {
381       CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data));
382     }
383   } else {
384     *(void **)data = NULL;
385   }
386   CeedCall(CeedQFunctionContextDestroy(&ctx));
387   return CEED_ERROR_SUCCESS;
388 }
389 
390 /**
391   @brief Restore context data of a `CeedQFunction`
392 
393   @param[in]     qf   `CeedQFunction`
394   @param[in,out] data Data to restore
395 
396   @return An error code: 0 - success, otherwise - failure
397 
398   @ref Backend
399 **/
400 int CeedQFunctionRestoreContextData(CeedQFunction qf, void *data) {
401   bool                 is_writable;
402   CeedQFunctionContext ctx;
403 
404   CeedCall(CeedQFunctionGetContext(qf, &ctx));
405   if (ctx) {
406     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
407     if (is_writable) {
408       CeedCall(CeedQFunctionContextRestoreData(ctx, data));
409     } else {
410       CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data));
411     }
412   }
413   CeedCall(CeedQFunctionContextDestroy(&ctx));
414   return CEED_ERROR_SUCCESS;
415 }
416 
417 /**
418   @brief Get true user context for a `CeedQFunction`
419 
420   Note: For all `CeedQFunction` this function will return the user `CeedQFunctionContext` and not interface context `CeedQFunctionContext`, if any such object exists.
421 
422   @param[in]  qf  `CeedQFunction`
423   @param[out] ctx Variable to store `CeedQFunctionContext`
424 
425   @return An error code: 0 - success, otherwise - failure
426   @ref Backend
427 **/
428 int CeedQFunctionGetInnerContext(CeedQFunction qf, CeedQFunctionContext *ctx) {
429   CeedQFunctionContext qf_ctx;
430 
431   CeedCall(CeedQFunctionGetContext(qf, &qf_ctx));
432   if (qf->is_fortran) {
433     CeedFortranContext fortran_ctx = NULL;
434 
435     CeedCall(CeedQFunctionContextGetData(qf_ctx, CEED_MEM_HOST, &fortran_ctx));
436     *ctx = fortran_ctx->inner_ctx;
437     CeedCall(CeedQFunctionContextRestoreData(qf_ctx, &fortran_ctx));
438   } else {
439     *ctx = qf_ctx;
440   }
441   CeedCall(CeedQFunctionContextDestroy(&qf_ctx));
442   return CEED_ERROR_SUCCESS;
443 }
444 
445 /**
446   @brief Get inner context data of a `CeedQFunction`
447 
448   @param[in]  qf       `CeedQFunction`
449   @param[in]  mem_type Memory type on which to access the data.
450                          If the backend uses a different memory type, this will perform a copy.
451   @param[out] data     Data on memory type mem_type
452 
453   @return An error code: 0 - success, otherwise - failure
454 
455   @ref Backend
456 **/
457 int CeedQFunctionGetInnerContextData(CeedQFunction qf, CeedMemType mem_type, void *data) {
458   bool                 is_writable;
459   CeedQFunctionContext ctx;
460 
461   CeedCall(CeedQFunctionGetInnerContext(qf, &ctx));
462   if (ctx) {
463     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
464     if (is_writable) {
465       CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data));
466     } else {
467       CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data));
468     }
469   } else {
470     *(void **)data = NULL;
471   }
472   return CEED_ERROR_SUCCESS;
473 }
474 
475 /**
476   @brief Restore inner context data of a `CeedQFunction`
477 
478   @param[in]     qf   `CeedQFunction`
479   @param[in,out] data Data to restore
480 
481   @return An error code: 0 - success, otherwise - failure
482 
483   @ref Backend
484 **/
485 int CeedQFunctionRestoreInnerContextData(CeedQFunction qf, void *data) {
486   bool                 is_writable;
487   CeedQFunctionContext ctx;
488 
489   CeedCall(CeedQFunctionGetInnerContext(qf, &ctx));
490   if (ctx) {
491     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
492     if (is_writable) {
493       CeedCall(CeedQFunctionContextRestoreData(ctx, data));
494     } else {
495       CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data));
496     }
497   }
498   *(void **)data = NULL;
499   return CEED_ERROR_SUCCESS;
500 }
501 
502 /**
503   @brief Determine if `CeedQFunction` is identity
504 
505   @param[in]  qf          `CeedQFunction`
506   @param[out] is_identity Variable to store identity status
507 
508   @return An error code: 0 - success, otherwise - failure
509 
510   @ref Backend
511 **/
512 int CeedQFunctionIsIdentity(CeedQFunction qf, bool *is_identity) {
513   *is_identity = qf->is_identity;
514   return CEED_ERROR_SUCCESS;
515 }
516 
517 /**
518   @brief Determine if `CeedQFunctionContext` is writable
519 
520   @param[in]  qf          `CeedQFunction`
521   @param[out] is_writable Variable to store context writeable status
522 
523   @return An error code: 0 - success, otherwise - failure
524 
525   @ref Backend
526 **/
527 int CeedQFunctionIsContextWritable(CeedQFunction qf, bool *is_writable) {
528   *is_writable = qf->is_context_writable;
529   return CEED_ERROR_SUCCESS;
530 }
531 
532 /**
533   @brief Get backend data of a `CeedQFunction`
534 
535   @param[in]  qf   `CeedQFunction`
536   @param[out] data Variable to store data
537 
538   @return An error code: 0 - success, otherwise - failure
539 
540   @ref Backend
541 **/
542 int CeedQFunctionGetData(CeedQFunction qf, void *data) {
543   *(void **)data = qf->data;
544   return CEED_ERROR_SUCCESS;
545 }
546 
547 /**
548   @brief Set backend data of a `CeedQFunction`
549 
550   @param[in,out] qf   `CeedQFunction`
551   @param[in]     data Data to set
552 
553   @return An error code: 0 - success, otherwise - failure
554 
555   @ref Backend
556 **/
557 int CeedQFunctionSetData(CeedQFunction qf, void *data) {
558   qf->data = data;
559   return CEED_ERROR_SUCCESS;
560 }
561 
562 /**
563   @brief Get a boolean value indicating if the `CeedQFunction` is immutable
564 
565   @param[in]  qf           `CeedOperator`
566   @param[out] is_immutable Variable to store immutability status
567 
568   @return An error code: 0 - success, otherwise - failure
569 
570   @ref Backend
571 **/
572 int CeedQFunctionIsImmutable(CeedQFunction qf, bool *is_immutable) {
573   *is_immutable = qf->is_immutable;
574   return CEED_ERROR_SUCCESS;
575 }
576 
577 /**
578   @brief Set the immutable flag of a `CeedQFunction` to `true`
579 
580   @param[in,out] qf `CeedQFunction`
581 
582   @return An error code: 0 - success, otherwise - failure
583 
584   @ref Backend
585 **/
586 int CeedQFunctionSetImmutable(CeedQFunction qf) {
587   qf->is_immutable = true;
588   return CEED_ERROR_SUCCESS;
589 }
590 
591 /**
592   @brief Increment the reference counter for a `CeedQFunction`
593 
594   @param[in,out] qf `CeedQFunction` to increment the reference counter
595 
596   @return An error code: 0 - success, otherwise - failure
597 
598   @ref Backend
599 **/
600 int CeedQFunctionReference(CeedQFunction qf) {
601   qf->ref_count++;
602   return CEED_ERROR_SUCCESS;
603 }
604 
605 /**
606   @brief Estimate number of FLOPs per quadrature required to apply `CeedQFunction`
607 
608   @param[in]  qf    `CeedQFunction` to estimate FLOPs for
609   @param[out] flops Address of variable to hold FLOPs estimate
610 
611   @ref Backend
612 **/
613 int CeedQFunctionGetFlopsEstimate(CeedQFunction qf, CeedSize *flops) {
614   *flops = qf->user_flop_estimate;
615   return CEED_ERROR_SUCCESS;
616 }
617 
618 /// @}
619 
620 /// ----------------------------------------------------------------------------
621 /// CeedQFunction Public API
622 /// ----------------------------------------------------------------------------
623 /// @addtogroup CeedQFunctionUser
624 /// @{
625 
626 /**
627   @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms
628 
629   @param[in]  ceed       `Ceed` object used to create the `CeedQFunction`
630   @param[in]  vec_length Vector length.
631                            Caller must ensure that number of quadrature points is a multiple of `vec_length`.
632   @param[in]  f          Function pointer to evaluate action at quadrature points.
633                            See `CeedQFunctionUser`.
634   @param[in]  source     Absolute path to source of `CeedQFunctionUser`, "\abs_path\file.h:function_name".
635                            The entire source file must only contain constructs supported by all targeted backends (i.e. CUDA for `/gpu/cuda`, OpenCL/SYCL for `/gpu/sycl`, etc.).
636                            The entire contents of this file and all locally included files are used during JiT compilation for GPU backends.
637                            The header `ceed/types.h` is preferred over `ceed.h` or `ceed/ceed.h` for `CeedQFunction` source files.
638                            The macro `CEED_RUNNING_JIT_PASS` is set during JiT and can be used to guard include statements that JiT compilers cannot use, such as `math.h` or `std*.h`.
639                            All source files must be at the provided filepath at runtime for JiT to function.
640   @param[out] qf         Address of the variable where the newly created `CeedQFunction` will be stored
641 
642   @return An error code: 0 - success, otherwise - failure
643 
644   See \ref CeedQFunctionUser for details on the call-back function `f` arguments.
645 
646   @ref User
647 **/
648 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vec_length, CeedQFunctionUser f, const char *source, CeedQFunction *qf) {
649   char *user_source_copy;
650 
651   if (!ceed->QFunctionCreate) {
652     Ceed delegate;
653 
654     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "QFunction"));
655     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement CeedQFunctionCreateInterior");
656     CeedCall(CeedQFunctionCreateInterior(delegate, vec_length, f, source, qf));
657     CeedCall(CeedDestroy(&delegate));
658     return CEED_ERROR_SUCCESS;
659   }
660 
661   CeedCheck(!strlen(source) || strrchr(source, ':'), ceed, CEED_ERROR_INCOMPLETE,
662             "Provided path to source does not include function name. Provided: \"%s\"\nRequired: \"\\abs_path\\file.h:function_name\"", source);
663 
664   CeedCall(CeedCalloc(1, qf));
665   CeedCall(CeedReferenceCopy(ceed, &(*qf)->ceed));
666   (*qf)->ref_count           = 1;
667   (*qf)->vec_length          = vec_length;
668   (*qf)->is_identity         = false;
669   (*qf)->is_context_writable = true;
670   (*qf)->function            = f;
671   (*qf)->user_flop_estimate  = -1;
672   if (strlen(source)) {
673     size_t user_source_len = strlen(source);
674 
675     CeedCall(CeedCalloc(user_source_len + 1, &user_source_copy));
676     memcpy(user_source_copy, source, user_source_len);
677     (*qf)->user_source = user_source_copy;
678   }
679   CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->input_fields));
680   CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->output_fields));
681   CeedCall(ceed->QFunctionCreate(*qf));
682   return CEED_ERROR_SUCCESS;
683 }
684 
685 /**
686   @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms by name
687 
688   @param[in]  ceed `Ceed` object used to create the `CeedQFunction`
689   @param[in]  name Name of `CeedQFunction` to use from gallery
690   @param[out] qf   Address of the variable where the newly created `CeedQFunction` will be stored
691 
692   @return An error code: 0 - success, otherwise - failure
693 
694   @ref User
695 **/
696 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, CeedQFunction *qf) {
697   size_t match_len = 0, match_index = UINT_MAX;
698 
699   CeedCall(CeedQFunctionRegisterAll());
700   // Find matching backend
701   CeedCheck(name, ceed, CEED_ERROR_INCOMPLETE, "No CeedQFunction name provided");
702   for (size_t i = 0; i < num_qfunctions; i++) {
703     size_t      n;
704     const char *curr_name = gallery_qfunctions[i].name;
705     for (n = 0; curr_name[n] && curr_name[n] == name[n]; n++) {
706     }
707     if (n > match_len) {
708       match_len   = n;
709       match_index = i;
710     }
711   }
712   CeedCheck(match_len > 0, ceed, CEED_ERROR_UNSUPPORTED, "No suitable gallery CeedQFunction");
713 
714   // Create QFunction
715   CeedCall(CeedQFunctionCreateInterior(ceed, gallery_qfunctions[match_index].vec_length, gallery_qfunctions[match_index].f,
716                                        gallery_qfunctions[match_index].source, qf));
717 
718   // QFunction specific setup
719   CeedCall(gallery_qfunctions[match_index].init(ceed, name, *qf));
720 
721   // Copy name
722   CeedCall(CeedStringAllocCopy(name, (char **)&(*qf)->gallery_name));
723   (*qf)->is_gallery = true;
724   return CEED_ERROR_SUCCESS;
725 }
726 
727 /**
728   @brief Create an identity `CeedQFunction`.
729 
730   Inputs are written into outputs in the order given.
731   This is useful for `CeedOperator that can be represented with only the action of a `CeedElemRestriction` and `CeedBasis`, such as restriction and prolongation operators for p-multigrid.
732   Backends may optimize `CeedOperator` with this `CeedQFunction` to avoid the copy of input data to output fields by using the same memory location for both.
733 
734   @param[in]  ceed     `Ceed` object used to create the `CeedQFunction`
735   @param[in]  size     Size of the `CeedQFunction` fields
736   @param[in]  in_mode  @ref CeedEvalMode for input to `CeedQFunction`
737   @param[in]  out_mode @ref CeedEvalMode for output to `CeedQFunction`
738   @param[out] qf       Address of the variable where the newly created `CeedQFunction` will be stored
739 
740   @return An error code: 0 - success, otherwise - failure
741 
742   @ref User
743 **/
744 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode in_mode, CeedEvalMode out_mode, CeedQFunction *qf) {
745   CeedQFunctionContext  ctx;
746   CeedContextFieldLabel size_label;
747 
748   CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Identity", qf));
749   CeedCall(CeedQFunctionAddInput(*qf, "input", size, in_mode));
750   CeedCall(CeedQFunctionAddOutput(*qf, "output", size, out_mode));
751 
752   (*qf)->is_identity = true;
753 
754   CeedCall(CeedQFunctionGetContext(*qf, &ctx));
755   CeedCall(CeedQFunctionContextGetFieldLabel(ctx, "size", &size_label));
756   CeedCall(CeedQFunctionContextSetInt32(ctx, size_label, &size));
757   CeedCall(CeedQFunctionContextDestroy(&ctx));
758   return CEED_ERROR_SUCCESS;
759 }
760 
761 /**
762   @brief Copy the pointer to a `CeedQFunction`.
763 
764   Both pointers should be destroyed with @ref CeedQFunctionDestroy().
765 
766   Note: If the value of `*qf_copy` passed to this function is non-NULL, then it is assumed that `*qf_copy` is a pointer to a `CeedQFunction`.
767         This `CeedQFunction` will be destroyed if `*qf_copy` is the only reference to this `CeedQFunction`.
768 
769   @param[in]  qf      `CeedQFunction` to copy reference to
770   @param[out] qf_copy Variable to store copied reference
771 
772   @return An error code: 0 - success, otherwise - failure
773 
774   @ref User
775 **/
776 int CeedQFunctionReferenceCopy(CeedQFunction qf, CeedQFunction *qf_copy) {
777   CeedCall(CeedQFunctionReference(qf));
778   CeedCall(CeedQFunctionDestroy(qf_copy));
779   *qf_copy = qf;
780   return CEED_ERROR_SUCCESS;
781 }
782 
783 /**
784   @brief Add a `CeedQFunction` input
785 
786   @param[in,out] qf         `CeedQFunction`
787   @param[in]     field_name Name of `CeedQFunction` field
788   @param[in]     size       Size of `CeedQFunction` field,
789                               (`num_comp * 1`) for @ref CEED_EVAL_NONE,
790                               (`num_comp * 1`) for @ref CEED_EVAL_INTERP for an \f$H^1\f$ space or (`num_comp * dim`) for an \f$H(\mathrm{div})\f$ or \f$H(\mathrm{curl})\f$ space,
791                               (`num_comp * dim`) for @ref CEED_EVAL_GRAD,
792                               (`num_comp * 1`) for @ref CEED_EVAL_DIV, and
793                               (`num_comp * curl_dim`) with `curl_dim = 1` if `dim < 3` otherwise `curl_dim = dim` for @ref CEED_EVAL_CURL.
794   @param[in]     eval_mode  @ref CEED_EVAL_NONE to use values directly,
795                               @ref CEED_EVAL_INTERP to use interpolated values,
796                               @ref CEED_EVAL_GRAD to use gradients,
797                               @ref CEED_EVAL_DIV to use divergence,
798                               @ref CEED_EVAL_CURL to use curl
799 
800   Note: In the user `CeedQFunctionUser`, the `in` argument list the fields in the order given by the calls to `CeedQFunctionAddInput`.
801 
802   @return An error code: 0 - success, otherwise - failure
803 
804   @ref User
805 **/
806 int CeedQFunctionAddInput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
807   bool is_immutable;
808 
809   CeedCall(CeedQFunctionIsImmutable(qf, &is_immutable));
810   CeedCheck(!is_immutable, CeedQFunctionReturnCeed(qf), CEED_ERROR_MAJOR, "QFunction cannot be changed after set as immutable");
811   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || size == 1, CeedQFunctionReturnCeed(qf), CEED_ERROR_DIMENSION, "CEED_EVAL_WEIGHT should have size 1");
812   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
813     CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR,
814               "CeedQFunction field names must be unique. Duplicate name: %s", field_name);
815   }
816   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
817     CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR,
818               "CeedQFunction field names must be unique. Duplicate name: %s", field_name);
819   }
820   CeedCall(CeedQFunctionFieldSet(&qf->input_fields[qf->num_input_fields], field_name, size, eval_mode));
821   qf->num_input_fields++;
822   return CEED_ERROR_SUCCESS;
823 }
824 
825 /**
826   @brief Add a `CeedQFunction` output
827 
828   @param[in,out] qf         `CeedQFunction`
829   @param[in]     field_name Name of `CeedQFunction` field
830   @param[in]     size       Size of `CeedQFunction` field,
831                               (`num_comp * 1`) for @ref CEED_EVAL_NONE,
832                               (`num_comp * 1`) for @ref CEED_EVAL_INTERP for an \f$H^1\f$ space or (`num_comp * dim`) for an \f$H(\mathrm{div})\f$ or \f$H(\mathrm{curl})\f$ space,
833                               (`num_comp * dim`) for @ref CEED_EVAL_GRAD,
834                               (`num_comp * 1`) for @ref CEED_EVAL_DIV, and
835                               (`num_comp * curl_dim`) with `curl_dim = 1` if `dim < 3` otherwise `curl_dim = dim` for @ref CEED_EVAL_CURL.
836   @param[in]     eval_mode  @ref CEED_EVAL_NONE to use values directly,
837                               @ref CEED_EVAL_INTERP to use interpolated values,
838                               @ref CEED_EVAL_GRAD to use gradients,
839                               @ref CEED_EVAL_DIV to use divergence,
840                               @ref CEED_EVAL_CURL to use curl.
841 
842   Note: In the user `CeedQFunctionUser`, the `out` argument list the fields in the order given by the calls to `CeedQFunctionAddOutput`.
843 
844   @return An error code: 0 - success, otherwise - failure
845 
846   @ref User
847 **/
848 int CeedQFunctionAddOutput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
849   bool is_immutable;
850 
851   CeedCall(CeedQFunctionIsImmutable(qf, &is_immutable));
852   CeedCheck(!is_immutable, CeedQFunctionReturnCeed(qf), CEED_ERROR_MAJOR, "CeedQFunction cannot be changed after set as immutable");
853   CeedCheck(eval_mode != CEED_EVAL_WEIGHT, CeedQFunctionReturnCeed(qf), CEED_ERROR_DIMENSION,
854             "Cannot create CeedQFunction output with CEED_EVAL_WEIGHT");
855   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
856     CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR,
857               "CeedQFunction field names must be unique");
858   }
859   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
860     CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR,
861               "CeedQFunction field names must be unique");
862   }
863   CeedCall(CeedQFunctionFieldSet(&qf->output_fields[qf->num_output_fields], field_name, size, eval_mode));
864   qf->num_output_fields++;
865   return CEED_ERROR_SUCCESS;
866 }
867 
868 /**
869   @brief Get the `CeedQFunctionField` of a `CeedQFunction`
870 
871   Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable.
872 
873   @param[in]  qf                `CeedQFunction`
874   @param[out] num_input_fields  Variable to store number of input fields
875   @param[out] input_fields      Variable to store input fields
876   @param[out] num_output_fields Variable to store number of output fields
877   @param[out] output_fields     Variable to store output fields
878 
879   @return An error code: 0 - success, otherwise - failure
880 
881   @ref Advanced
882 **/
883 int CeedQFunctionGetFields(CeedQFunction qf, CeedInt *num_input_fields, CeedQFunctionField **input_fields, CeedInt *num_output_fields,
884                            CeedQFunctionField **output_fields) {
885   CeedCall(CeedQFunctionSetImmutable(qf));
886   if (num_input_fields) *num_input_fields = qf->num_input_fields;
887   if (input_fields) *input_fields = qf->input_fields;
888   if (num_output_fields) *num_output_fields = qf->num_output_fields;
889   if (output_fields) *output_fields = qf->output_fields;
890   return CEED_ERROR_SUCCESS;
891 }
892 
893 /**
894   @brief Get the name of a `CeedQFunctionField`
895 
896   @param[in]  qf_field   `CeedQFunctionField`
897   @param[out] field_name Variable to store the field name
898 
899   @return An error code: 0 - success, otherwise - failure
900 
901   @ref Advanced
902 **/
903 int CeedQFunctionFieldGetName(CeedQFunctionField qf_field, const char **field_name) {
904   *field_name = qf_field->field_name;
905   return CEED_ERROR_SUCCESS;
906 }
907 
908 /**
909   @brief Get the number of components of a `CeedQFunctionField`
910 
911   @param[in]  qf_field `CeedQFunctionField`
912   @param[out] size     Variable to store the size of the field
913 
914   @return An error code: 0 - success, otherwise - failure
915 
916   @ref Advanced
917 **/
918 int CeedQFunctionFieldGetSize(CeedQFunctionField qf_field, CeedInt *size) {
919   *size = qf_field->size;
920   return CEED_ERROR_SUCCESS;
921 }
922 
923 /**
924   @brief Get the @ref CeedEvalMode of a `CeedQFunctionField`
925 
926   @param[in]  qf_field  `CeedQFunctionField`
927   @param[out] eval_mode Variable to store the field evaluation mode
928 
929   @return An error code: 0 - success, otherwise - failure
930 
931   @ref Advanced
932 **/
933 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qf_field, CeedEvalMode *eval_mode) {
934   *eval_mode = qf_field->eval_mode;
935   return CEED_ERROR_SUCCESS;
936 }
937 
938 /**
939   @brief Get the data of a `CeedQFunctionField`.
940 
941   Any arguments set as `NULL` are ignored.
942 
943   @param[in]  qf_field   `CeedQFunctionField`
944   @param[out] field_name Variable to store the field name
945   @param[out] size       Variable to store the size of the field
946   @param[out] eval_mode  Variable to store the field evaluation mode
947 
948   @return An error code: 0 - success, otherwise - failure
949 
950   @ref Advanced
951 **/
952 int CeedQFunctionFieldGetData(CeedQFunctionField qf_field, const char **field_name, CeedInt *size, CeedEvalMode *eval_mode) {
953   if (field_name) CeedCall(CeedQFunctionFieldGetName(qf_field, field_name));
954   if (size) CeedCall(CeedQFunctionFieldGetSize(qf_field, size));
955   if (eval_mode) CeedCall(CeedQFunctionFieldGetEvalMode(qf_field, eval_mode));
956   return CEED_ERROR_SUCCESS;
957 }
958 
959 /**
960   @brief Set global context for a `CeedQFunction`
961 
962   @param[in,out] qf  `CeedQFunction`
963   @param[in]     ctx Context data to set
964 
965   @return An error code: 0 - success, otherwise - failure
966 
967   @ref User
968 **/
969 int CeedQFunctionSetContext(CeedQFunction qf, CeedQFunctionContext ctx) {
970   CeedCall(CeedQFunctionContextDestroy(&qf->ctx));
971   qf->ctx = ctx;
972   if (ctx) CeedCall(CeedQFunctionContextReference(ctx));
973   return CEED_ERROR_SUCCESS;
974 }
975 
976 /**
977   @brief Set writability of `CeedQFunctionContext` when calling the `CeedQFunctionUser`.
978 
979   The default value is `is_writable == true`.
980 
981   Setting `is_writable == true` indicates the `CeedQFunctionUser` writes into the `CeedQFunctionContext` and requires memory synchronization after calling @ref CeedQFunctionApply().
982 
983   Setting 'is_writable == false' asserts that `CeedQFunctionUser` does not modify the `CeedQFunctionContext`.
984   Violating this assertion may lead to inconsistent data.
985 
986   Setting `is_writable == false` may offer a performance improvement on GPU backends.
987 
988   @param[in,out] qf          `CeedQFunction`
989   @param[in]     is_writable Boolean flag for writability status
990 
991   @return An error code: 0 - success, otherwise - failure
992 
993   @ref User
994 **/
995 int CeedQFunctionSetContextWritable(CeedQFunction qf, bool is_writable) {
996   qf->is_context_writable = is_writable;
997   return CEED_ERROR_SUCCESS;
998 }
999 
1000 /**
1001   @brief Set estimated number of FLOPs per quadrature required to apply `CeedQFunction`
1002 
1003   @param[in]  qf    `CeedQFunction` to estimate FLOPs for
1004   @param[out] flops FLOPs per quadrature point estimate
1005 
1006   @ref Backend
1007 **/
1008 int CeedQFunctionSetUserFlopsEstimate(CeedQFunction qf, CeedSize flops) {
1009   CeedCheck(flops >= 0, CeedQFunctionReturnCeed(qf), CEED_ERROR_INCOMPATIBLE, "Must set non-negative FLOPs estimate");
1010   qf->user_flop_estimate = flops;
1011   return CEED_ERROR_SUCCESS;
1012 }
1013 
1014 /**
1015   @brief Set the number of tabs to indent for @ref CeedQFunctionView() output
1016 
1017   @param[in] qf       `CeedQFunction` to set the number of view tabs
1018   @param[in] num_tabs Number of view tabs to set
1019 
1020   @return Error code: 0 - success, otherwise - failure
1021 
1022   @ref User
1023 **/
1024 int CeedQFunctionSetNumViewTabs(CeedQFunction qf, CeedInt num_tabs) {
1025   CeedCheck(num_tabs >= 0, CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR, "Number of view tabs must be non-negative");
1026   qf->num_tabs = num_tabs;
1027   return CEED_ERROR_SUCCESS;
1028 }
1029 
1030 /**
1031   @brief Get the number of tabs to indent for @ref CeedQFunctionView() output
1032 
1033   @param[in]  qf       `CeedQFunction` to get the number of view tabs
1034   @param[out] num_tabs Number of view tabs
1035 
1036   @return Error code: 0 - success, otherwise - failure
1037 
1038   @ref User
1039 **/
1040 int CeedQFunctionGetNumViewTabs(CeedQFunction qf, CeedInt *num_tabs) {
1041   *num_tabs = qf->num_tabs;
1042   return CEED_ERROR_SUCCESS;
1043 }
1044 
1045 /**
1046   @brief View a `CeedQFunction`
1047 
1048   @param[in] qf     `CeedQFunction` to view
1049   @param[in] stream Stream to write; typically `stdout` or a file
1050 
1051   @return Error code: 0 - success, otherwise - failure
1052 
1053   @ref User
1054 **/
1055 int CeedQFunctionView(CeedQFunction qf, FILE *stream) {
1056   char       *tabs = NULL;
1057   const char *name;
1058 
1059   {
1060     CeedInt num_tabs = 0;
1061 
1062     CeedCall(CeedQFunctionGetNumViewTabs(qf, &num_tabs));
1063     CeedCall(CeedCalloc(CEED_TAB_WIDTH * num_tabs + 1, &tabs));
1064     for (CeedInt i = 0; i < CEED_TAB_WIDTH * num_tabs; i++) tabs[i] = ' ';
1065   }
1066 
1067   CeedCall(CeedQFunctionGetName(qf, &name));
1068   fprintf(stream, "%s%sCeedQFunction - %s\n", tabs, qf->is_gallery ? "Gallery " : "User ", name);
1069 
1070   fprintf(stream, "%s  %" CeedInt_FMT " input field%s:\n", tabs, qf->num_input_fields, qf->num_input_fields > 1 ? "s" : "");
1071   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
1072     CeedCall(CeedQFunctionFieldView(qf->input_fields[i], i, 1, tabs, stream));
1073   }
1074 
1075   fprintf(stream, "%s  %" CeedInt_FMT " output field%s:\n", tabs, qf->num_output_fields, qf->num_output_fields > 1 ? "s" : "");
1076   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
1077     CeedCall(CeedQFunctionFieldView(qf->output_fields[i], i, 0, tabs, stream));
1078   }
1079   CeedCall(CeedFree(&tabs));
1080   return CEED_ERROR_SUCCESS;
1081 }
1082 
1083 /**
1084   @brief Get the `Ceed` associated with a `CeedQFunction`
1085 
1086   @param[in]  qf   `CeedQFunction`
1087   @param[out] ceed Variable to store`Ceed`
1088 
1089   @return An error code: 0 - success, otherwise - failure
1090 
1091   @ref Advanced
1092 **/
1093 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) {
1094   *ceed = NULL;
1095   CeedCall(CeedReferenceCopy(CeedQFunctionReturnCeed(qf), ceed));
1096   return CEED_ERROR_SUCCESS;
1097 }
1098 
1099 /**
1100   @brief Return the `Ceed` associated with a `CeedQFunction`
1101 
1102   @param[in]  qf   `CeedQFunction`
1103 
1104   @return `Ceed` associated with the `qf`
1105 
1106   @ref Advanced
1107 **/
1108 Ceed CeedQFunctionReturnCeed(CeedQFunction qf) { return qf->ceed; }
1109 
1110 /**
1111   @brief Apply the action of a `CeedQFunction`
1112 
1113   Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable.
1114 
1115   @param[in]  qf `CeedQFunction`
1116   @param[in]  Q  Number of quadrature points
1117   @param[in]  u  Array of input `CeedVector`
1118   @param[out] v  Array of output `CeedVector`
1119 
1120   @return An error code: 0 - success, otherwise - failure
1121 
1122   @ref User
1123 **/
1124 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, CeedVector *u, CeedVector *v) {
1125   CeedInt vec_length;
1126 
1127   CeedCheck(qf->Apply, CeedQFunctionReturnCeed(qf), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedQFunctionApply");
1128   CeedCall(CeedQFunctionGetVectorLength(qf, &vec_length));
1129   CeedCheck(Q % vec_length == 0, CeedQFunctionReturnCeed(qf), CEED_ERROR_DIMENSION,
1130             "Number of quadrature points %" CeedInt_FMT " must be a multiple of %" CeedInt_FMT, Q, qf->vec_length);
1131   CeedCall(CeedQFunctionSetImmutable(qf));
1132   CeedCall(qf->Apply(qf, Q, u, v));
1133   return CEED_ERROR_SUCCESS;
1134 }
1135 
1136 /**
1137   @brief Destroy a `CeedQFunction`
1138 
1139   @param[in,out] qf `CeedQFunction` to destroy
1140 
1141   @return An error code: 0 - success, otherwise - failure
1142 
1143   @ref User
1144 **/
1145 int CeedQFunctionDestroy(CeedQFunction *qf) {
1146   if (!*qf || --(*qf)->ref_count > 0) {
1147     *qf = NULL;
1148     return CEED_ERROR_SUCCESS;
1149   }
1150   // Backend destroy
1151   if ((*qf)->Destroy) {
1152     CeedCall((*qf)->Destroy(*qf));
1153   }
1154   // Free fields
1155   for (CeedInt i = 0; i < (*qf)->num_input_fields; i++) {
1156     CeedCall(CeedFree(&(*(*qf)->input_fields[i]).field_name));
1157     CeedCall(CeedFree(&(*qf)->input_fields[i]));
1158   }
1159   for (CeedInt i = 0; i < (*qf)->num_output_fields; i++) {
1160     CeedCall(CeedFree(&(*(*qf)->output_fields[i]).field_name));
1161     CeedCall(CeedFree(&(*qf)->output_fields[i]));
1162   }
1163   CeedCall(CeedFree(&(*qf)->input_fields));
1164   CeedCall(CeedFree(&(*qf)->output_fields));
1165 
1166   // User context data object
1167   CeedCall(CeedQFunctionContextDestroy(&(*qf)->ctx));
1168 
1169   CeedCall(CeedFree(&(*qf)->user_source));
1170   CeedCall(CeedFree(&(*qf)->source_path));
1171   CeedCall(CeedFree(&(*qf)->gallery_name));
1172   CeedCall(CeedFree(&(*qf)->kernel_name));
1173   CeedCall(CeedDestroy(&(*qf)->ceed));
1174   CeedCall(CeedFree(qf));
1175   return CEED_ERROR_SUCCESS;
1176 }
1177 
1178 /// @}
1179