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