xref: /libCEED/interface/ceed-qfunction.c (revision ca7355308a14805110a7d98dabf1f658d5ec16d5)
1 // Copyright (c) 2017-2022, 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 QFunction 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 QFunction
51 
52   @param[in]  name       Name for this backend to respond to
53   @param[in]  source     Absolute path to source of QFunction, "\path\CEED_DIR\gallery\folder\file.h:function_name"
54   @param[in]  vec_length Vector length. Caller must ensure that number of quadrature points is a multiple of vec_length.
55   @param[in]  f          Function pointer to evaluate action at quadrature points.
56                            See \ref CeedQFunctionUser.
57   @param[in]  init       Initialization function called by CeedQFunctionInit() when the QFunction is selected.
58 
59   @return An error code: 0 - success, otherwise - failure
60 
61   @ref Developer
62 **/
63 int CeedQFunctionRegister(const char *name, const char *source, CeedInt vec_length, CeedQFunctionUser f,
64                           int (*init)(Ceed, const char *, CeedQFunction)) {
65   const char *relative_file_path;
66 
67   CeedCheck(num_qfunctions < sizeof(gallery_qfunctions) / sizeof(gallery_qfunctions[0]), NULL, CEED_ERROR_MAJOR, "Too many gallery QFunctions");
68 
69   CeedDebugEnv("Gallery Register: %s", name);
70 
71   CeedCall(CeedGetJitRelativePath(source, &relative_file_path));
72 
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   return CEED_ERROR_SUCCESS;
82 }
83 
84 /**
85   @brief Set a CeedQFunction field, used by CeedQFunctionAddInput/Output
86 
87   @param[out] f           CeedQFunctionField
88   @param[in]  field_name  Name of QFunction field
89   @param[in]  size        Size of QFunction field, (num_comp * 1) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_WEIGHT,
90 (num_comp * 1) for @ref CEED_EVAL_INTERP for an H^1 space or (num_comp * dim) for an H(div) or H(curl) space,
91 (num_comp * dim) for @ref CEED_EVAL_GRAD, or (num_comp * 1) for @ref CEED_EVAL_DIV, and
92 (num_comp * curl_dim) with curl_dim = 1 if dim < 3 else dim for @ref CEED_EVAL_CURL.
93   @param[in]  eval_mode   \ref CEED_EVAL_NONE to use values directly,
94                             \ref CEED_EVAL_WEIGHT to use quadrature weights,
95                             \ref CEED_EVAL_INTERP to use interpolated values,
96                             \ref CEED_EVAL_GRAD to use gradients,
97                             \ref CEED_EVAL_DIV to use divergence,
98                             \ref CEED_EVAL_CURL to use curl.
99 
100   @return An error code: 0 - success, otherwise - failure
101 
102   @ref Developer
103 **/
104 static int CeedQFunctionFieldSet(CeedQFunctionField *f, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
105   CeedCall(CeedCalloc(1, f));
106   CeedCall(CeedStringAllocCopy(field_name, (char **)&(*f)->field_name));
107   (*f)->size      = size;
108   (*f)->eval_mode = eval_mode;
109   return CEED_ERROR_SUCCESS;
110 }
111 
112 /**
113   @brief View a field of a CeedQFunction
114 
115   @param[in] field        QFunction field to view
116   @param[in] field_number Number of field being viewed
117   @param[in] in           true for input field, false for output
118   @param[in] stream       Stream to view to, e.g., stdout
119 
120   @return An error code: 0 - success, otherwise - failure
121 
122   @ref Utility
123 **/
124 static int CeedQFunctionFieldView(CeedQFunctionField field, CeedInt field_number, bool in, FILE *stream) {
125   const char  *inout = in ? "Input" : "Output";
126   char        *field_name;
127   CeedInt      size;
128   CeedEvalMode eval_mode;
129 
130   CeedCall(CeedQFunctionFieldGetName(field, &field_name));
131   CeedCall(CeedQFunctionFieldGetSize(field, &size));
132   CeedCall(CeedQFunctionFieldGetEvalMode(field, &eval_mode));
133   fprintf(stream,
134           "    %s field %" CeedInt_FMT
135           ":\n"
136           "      Name: \"%s\"\n"
137           "      Size: %" CeedInt_FMT
138           "\n"
139           "      EvalMode: \"%s\"\n",
140           inout, field_number, field_name, size, CeedEvalModes[eval_mode]);
141   return CEED_ERROR_SUCCESS;
142 }
143 
144 /**
145   @brief Set flag to determine if Fortran interface is used
146 
147   @param[in,out] qf     CeedQFunction
148   @param[in]     status Boolean value to set as Fortran status
149 
150   @return An error code: 0 - success, otherwise - failure
151 
152   @ref Backend
153 **/
154 int CeedQFunctionSetFortranStatus(CeedQFunction qf, bool status) {
155   qf->is_fortran = status;
156   return CEED_ERROR_SUCCESS;
157 }
158 
159 /// @}
160 
161 /// ----------------------------------------------------------------------------
162 /// CeedQFunction Backend API
163 /// ----------------------------------------------------------------------------
164 /// @addtogroup CeedQFunctionBackend
165 /// @{
166 
167 /**
168   @brief Get the vector length of a CeedQFunction
169 
170   @param[in]  qf         CeedQFunction
171   @param[out] vec_length Variable to store vector length
172 
173   @return An error code: 0 - success, otherwise - failure
174 
175   @ref Backend
176 **/
177 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vec_length) {
178   *vec_length = qf->vec_length;
179   return CEED_ERROR_SUCCESS;
180 }
181 
182 /**
183   @brief Get the number of inputs and outputs to a CeedQFunction
184 
185   @param[in]  qf         CeedQFunction
186   @param[out] num_input  Variable to store number of input fields
187   @param[out] num_output Variable to store number of output fields
188 
189   @return An error code: 0 - success, otherwise - failure
190 
191   @ref Backend
192 **/
193 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *num_input, CeedInt *num_output) {
194   if (num_input) *num_input = qf->num_input_fields;
195   if (num_output) *num_output = qf->num_output_fields;
196   return CEED_ERROR_SUCCESS;
197 }
198 
199 /**
200   @brief Get the name of the user function for a CeedQFunction
201 
202   @param[in]  qf          CeedQFunction
203   @param[out] kernel_name Variable to store source path string
204 
205   @return An error code: 0 - success, otherwise - failure
206 
207   @ref Backend
208 **/
209 int CeedQFunctionGetKernelName(CeedQFunction qf, char **kernel_name) {
210   if (!qf->kernel_name) {
211     Ceed  ceed;
212     char *kernel_name_copy;
213 
214     CeedCall(CeedQFunctionGetCeed(qf, &ceed));
215     if (qf->user_source) {
216       const char *kernel_name     = strrchr(qf->user_source, ':') + 1;
217       size_t      kernel_name_len = strlen(kernel_name);
218 
219       CeedCall(CeedCalloc(kernel_name_len + 1, &kernel_name_copy));
220       memcpy(kernel_name_copy, kernel_name, kernel_name_len);
221     } else {
222       CeedCall(CeedCalloc(1, &kernel_name_copy));
223     }
224     qf->kernel_name = kernel_name_copy;
225   }
226 
227   *kernel_name = (char *)qf->kernel_name;
228   return CEED_ERROR_SUCCESS;
229 }
230 
231 /**
232   @brief Get the source path string for a CeedQFunction
233 
234   @param[in]  qf          CeedQFunction
235   @param[out] source_path Variable to store source path string
236 
237   @return An error code: 0 - success, otherwise - failure
238 
239   @ref Backend
240 **/
241 int CeedQFunctionGetSourcePath(CeedQFunction qf, char **source_path) {
242   if (!qf->source_path && qf->user_source) {
243     Ceed        ceed;
244     bool        is_absolute_path;
245     char       *absolute_path, *source_path_copy;
246     const char *kernel_name     = strrchr(qf->user_source, ':') + 1;
247     size_t      kernel_name_len = strlen(kernel_name);
248 
249     CeedCall(CeedQFunctionGetCeed(qf, &ceed));
250     CeedCall(CeedCheckFilePath(ceed, qf->user_source, &is_absolute_path));
251     if (is_absolute_path) {
252       absolute_path = (char *)qf->user_source;
253     } else {
254       CeedCall(CeedGetJitAbsolutePath(ceed, qf->user_source, &absolute_path));
255     }
256 
257     size_t source_len = strlen(absolute_path) - kernel_name_len - 1;
258 
259     CeedCall(CeedCalloc(source_len + 1, &source_path_copy));
260     memcpy(source_path_copy, absolute_path, source_len);
261     qf->source_path = source_path_copy;
262 
263     if (!is_absolute_path) CeedCall(CeedFree(&absolute_path));
264   }
265 
266   *source_path = (char *)qf->source_path;
267   return CEED_ERROR_SUCCESS;
268 }
269 
270 /**
271   @brief Initialize and load QFunction source file into string buffer, including full text of local files in place of `#include "local.h"`.
272 
273   The `buffer` is set to `NULL` if there is no QFunction source file.
274 
275   Note: Caller is responsible for freeing the string buffer with `CeedFree()`.
276 
277   @param[in]  qf            CeedQFunction
278   @param[out] source_buffer String buffer for source file contents
279 
280   @return An error code: 0 - success, otherwise - failure
281 
282   @ref Backend
283 **/
284 int CeedQFunctionLoadSourceToBuffer(CeedQFunction qf, char **source_buffer) {
285   char *source_path;
286 
287   CeedCall(CeedQFunctionGetSourcePath(qf, &source_path));
288   *source_buffer = NULL;
289   if (source_path) {
290     CeedCall(CeedLoadSourceToBuffer(qf->ceed, source_path, source_buffer));
291   }
292   return CEED_ERROR_SUCCESS;
293 }
294 
295 /**
296   @brief Get the User Function for a CeedQFunction
297 
298   @param[in]  qf CeedQFunction
299   @param[out] f  Variable to store user function
300 
301   @return An error code: 0 - success, otherwise - failure
302 
303   @ref Backend
304 **/
305 int CeedQFunctionGetUserFunction(CeedQFunction qf, CeedQFunctionUser *f) {
306   *f = qf->function;
307   return CEED_ERROR_SUCCESS;
308 }
309 
310 /**
311   @brief Get global context for a CeedQFunction.
312 
313   Note: For QFunctions from the Fortran interface, this function will return the Fortran context CeedQFunctionContext.
314 
315   @param[in]  qf  CeedQFunction
316   @param[out] ctx Variable to store CeedQFunctionContext
317 
318   @return An error code: 0 - success, otherwise - failure
319 
320   @ref Backend
321 **/
322 int CeedQFunctionGetContext(CeedQFunction qf, CeedQFunctionContext *ctx) {
323   *ctx = qf->ctx;
324   return CEED_ERROR_SUCCESS;
325 }
326 
327 /**
328   @brief Get context data of a CeedQFunction
329 
330   @param[in]  qf       CeedQFunction
331   @param[in]  mem_type Memory type on which to access the data.
332                          If the backend uses a different memory type, this will perform a copy.
333   @param[out] data     Data on memory type mem_type
334 
335   @return An error code: 0 - success, otherwise - failure
336 
337   @ref Backend
338 **/
339 int CeedQFunctionGetContextData(CeedQFunction qf, CeedMemType mem_type, void *data) {
340   bool                 is_writable;
341   CeedQFunctionContext ctx;
342 
343   CeedCall(CeedQFunctionGetContext(qf, &ctx));
344   if (ctx) {
345     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
346     if (is_writable) {
347       CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data));
348     } else {
349       CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data));
350     }
351   } else {
352     *(void **)data = NULL;
353   }
354   return CEED_ERROR_SUCCESS;
355 }
356 
357 /**
358   @brief Restore context data of a CeedQFunction
359 
360   @param[in]     qf   CeedQFunction
361   @param[in,out] data Data to restore
362 
363   @return An error code: 0 - success, otherwise - failure
364 
365   @ref Backend
366 **/
367 int CeedQFunctionRestoreContextData(CeedQFunction qf, void *data) {
368   bool                 is_writable;
369   CeedQFunctionContext ctx;
370 
371   CeedCall(CeedQFunctionGetContext(qf, &ctx));
372   if (ctx) {
373     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
374     if (is_writable) {
375       CeedCall(CeedQFunctionContextRestoreData(ctx, data));
376     } else {
377       CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data));
378     }
379   }
380   *(void **)data = NULL;
381   return CEED_ERROR_SUCCESS;
382 }
383 
384 /**
385   @brief Get true user context for a CeedQFunction
386 
387   Note: For all QFunctions this function will return the user CeedQFunctionContext and not interface context CeedQFunctionContext, if any
388 such object exists.
389 
390   @param[in]  qf  CeedQFunction
391   @param[out] ctx Variable to store CeedQFunctionContext
392 
393   @return An error code: 0 - success, otherwise - failure
394   @ref Backend
395 **/
396 int CeedQFunctionGetInnerContext(CeedQFunction qf, CeedQFunctionContext *ctx) {
397   if (qf->is_fortran) {
398     CeedFortranContext fortran_ctx = NULL;
399 
400     CeedCall(CeedQFunctionContextGetData(qf->ctx, CEED_MEM_HOST, &fortran_ctx));
401     *ctx = fortran_ctx->inner_ctx;
402     CeedCall(CeedQFunctionContextRestoreData(qf->ctx, (void *)&fortran_ctx));
403   } else {
404     *ctx = qf->ctx;
405   }
406   return CEED_ERROR_SUCCESS;
407 }
408 
409 /**
410   @brief Get inner context data of a CeedQFunction
411 
412   @param[in]  qf       CeedQFunction
413   @param[in]  mem_type Memory type on which to access the data.
414                          If the backend uses a different memory type, this will perform a copy.
415   @param[out] data     Data on memory type mem_type
416 
417   @return An error code: 0 - success, otherwise - failure
418 
419   @ref Backend
420 **/
421 int CeedQFunctionGetInnerContextData(CeedQFunction qf, CeedMemType mem_type, void *data) {
422   bool                 is_writable;
423   CeedQFunctionContext ctx;
424 
425   CeedCall(CeedQFunctionGetInnerContext(qf, &ctx));
426   if (ctx) {
427     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
428     if (is_writable) {
429       CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data));
430     } else {
431       CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data));
432     }
433   } else {
434     *(void **)data = NULL;
435   }
436   return CEED_ERROR_SUCCESS;
437 }
438 
439 /**
440   @brief Restore inner context data of a CeedQFunction
441 
442   @param[in]     qf   CeedQFunction
443   @param[in,out] data Data to restore
444 
445   @return An error code: 0 - success, otherwise - failure
446 
447   @ref Backend
448 **/
449 int CeedQFunctionRestoreInnerContextData(CeedQFunction qf, void *data) {
450   bool                 is_writable;
451   CeedQFunctionContext ctx;
452 
453   CeedCall(CeedQFunctionGetInnerContext(qf, &ctx));
454   if (ctx) {
455     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
456     if (is_writable) {
457       CeedCall(CeedQFunctionContextRestoreData(ctx, data));
458     } else {
459       CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data));
460     }
461   }
462   *(void **)data = NULL;
463   return CEED_ERROR_SUCCESS;
464 }
465 
466 /**
467   @brief Determine if QFunction is identity
468 
469   @param[in]  qf          CeedQFunction
470   @param[out] is_identity Variable to store identity status
471 
472   @return An error code: 0 - success, otherwise - failure
473 
474   @ref Backend
475 **/
476 int CeedQFunctionIsIdentity(CeedQFunction qf, bool *is_identity) {
477   *is_identity = qf->is_identity;
478   return CEED_ERROR_SUCCESS;
479 }
480 
481 /**
482   @brief Determine if QFunctionContext is writable
483 
484   @param[in]  qf          CeedQFunction
485   @param[out] is_writable Variable to store context writeable status
486 
487   @return An error code: 0 - success, otherwise - failure
488 
489   @ref Backend
490 **/
491 int CeedQFunctionIsContextWritable(CeedQFunction qf, bool *is_writable) {
492   *is_writable = qf->is_context_writable;
493   return CEED_ERROR_SUCCESS;
494 }
495 
496 /**
497   @brief Get backend data of a CeedQFunction
498 
499   @param[in]  qf   CeedQFunction
500   @param[out] data Variable to store data
501 
502   @return An error code: 0 - success, otherwise - failure
503 
504   @ref Backend
505 **/
506 int CeedQFunctionGetData(CeedQFunction qf, void *data) {
507   *(void **)data = qf->data;
508   return CEED_ERROR_SUCCESS;
509 }
510 
511 /**
512   @brief Set backend data of a CeedQFunction
513 
514   @param[in,out] qf   CeedQFunction
515   @param[in]     data Data to set
516 
517   @return An error code: 0 - success, otherwise - failure
518 
519   @ref Backend
520 **/
521 int CeedQFunctionSetData(CeedQFunction qf, void *data) {
522   qf->data = data;
523   return CEED_ERROR_SUCCESS;
524 }
525 
526 /**
527   @brief Increment the reference counter for a CeedQFunction
528 
529   @param[in,out] qf CeedQFunction to increment the reference counter
530 
531   @return An error code: 0 - success, otherwise - failure
532 
533   @ref Backend
534 **/
535 int CeedQFunctionReference(CeedQFunction qf) {
536   qf->ref_count++;
537   return CEED_ERROR_SUCCESS;
538 }
539 
540 /**
541   @brief Estimate number of FLOPs per quadrature required to apply QFunction
542 
543   @param[in]  qf    QFunction to estimate FLOPs for
544   @param[out] flops Address of variable to hold FLOPs estimate
545 
546   @ref Backend
547 **/
548 int CeedQFunctionGetFlopsEstimate(CeedQFunction qf, CeedSize *flops) {
549   CeedCheck(qf->user_flop_estimate > -1, qf->ceed, CEED_ERROR_INCOMPLETE, "Must set FLOPs estimate with CeedQFunctionSetUserFlopsEstimate");
550   *flops = qf->user_flop_estimate;
551   return CEED_ERROR_SUCCESS;
552 }
553 
554 /// @}
555 
556 /// ----------------------------------------------------------------------------
557 /// CeedQFunction Public API
558 /// ----------------------------------------------------------------------------
559 /// @addtogroup CeedQFunctionUser
560 /// @{
561 
562 /**
563   @brief Create a CeedQFunction for evaluating interior (volumetric) terms.
564 
565   @param[in]  ceed       Ceed object where the CeedQFunction will be created
566   @param[in]  vec_length Vector length. Caller must ensure that number of quadrature points is a multiple of vec_length.
567   @param[in]  f          Function pointer to evaluate action at quadrature points.
568                            See \ref CeedQFunctionUser.
569   @param[in]  source     Absolute path to source of QFunction, "\abs_path\file.h:function_name".
570                            For support across all backends, this source must only contain constructs supported by C99, C++11, and CUDA.
571   @param[out] qf         Address of the variable where the newly created CeedQFunction will be stored
572 
573   @return An error code: 0 - success, otherwise - failure
574 
575   See \ref CeedQFunctionUser for details on the call-back function @a f's arguments.
576 
577   @ref User
578 **/
579 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vec_length, CeedQFunctionUser f, const char *source, CeedQFunction *qf) {
580   char *user_source_copy;
581 
582   if (!ceed->QFunctionCreate) {
583     Ceed delegate;
584 
585     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "QFunction"));
586     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support QFunctionCreate");
587     CeedCall(CeedQFunctionCreateInterior(delegate, vec_length, f, source, qf));
588     return CEED_ERROR_SUCCESS;
589   }
590 
591   CeedCheck(!strlen(source) || strrchr(source, ':'), ceed, CEED_ERROR_INCOMPLETE,
592             "Provided path to source does not include function name. Provided: \"%s\"\nRequired: \"\\abs_path\\file.h:function_name\"", source);
593 
594   CeedCall(CeedCalloc(1, qf));
595   CeedCall(CeedReferenceCopy(ceed, &(*qf)->ceed));
596   (*qf)->ref_count           = 1;
597   (*qf)->vec_length          = vec_length;
598   (*qf)->is_identity         = false;
599   (*qf)->is_context_writable = true;
600   (*qf)->function            = f;
601   (*qf)->user_flop_estimate  = -1;
602   if (strlen(source)) {
603     size_t user_source_len = strlen(source);
604 
605     CeedCall(CeedCalloc(user_source_len + 1, &user_source_copy));
606     memcpy(user_source_copy, source, user_source_len);
607     (*qf)->user_source = user_source_copy;
608   }
609   CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->input_fields));
610   CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->output_fields));
611   CeedCall(ceed->QFunctionCreate(*qf));
612   return CEED_ERROR_SUCCESS;
613 }
614 
615 /**
616   @brief Create a CeedQFunction for evaluating interior (volumetric) terms by name.
617 
618   @param[in]  ceed Ceed object where the CeedQFunction will be created
619   @param[in]  name Name of QFunction to use from gallery
620   @param[out] qf   Address of the variable where the newly created CeedQFunction will be stored
621 
622   @return An error code: 0 - success, otherwise - failure
623 
624   @ref User
625 **/
626 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, CeedQFunction *qf) {
627   size_t match_len = 0, match_index = UINT_MAX;
628 
629   CeedCall(CeedQFunctionRegisterAll());
630   // Find matching backend
631   CeedCheck(name, ceed, CEED_ERROR_INCOMPLETE, "No QFunction name provided");
632   for (size_t i = 0; i < num_qfunctions; i++) {
633     size_t      n;
634     const char *curr_name = gallery_qfunctions[i].name;
635     for (n = 0; curr_name[n] && curr_name[n] == name[n]; n++) {
636     }
637     if (n > match_len) {
638       match_len   = n;
639       match_index = i;
640     }
641   }
642   CeedCheck(match_len > 0, ceed, CEED_ERROR_UNSUPPORTED, "No suitable gallery QFunction");
643 
644   // Create QFunction
645   CeedCall(CeedQFunctionCreateInterior(ceed, gallery_qfunctions[match_index].vec_length, gallery_qfunctions[match_index].f,
646                                        gallery_qfunctions[match_index].source, qf));
647 
648   // QFunction specific setup
649   CeedCall(gallery_qfunctions[match_index].init(ceed, name, *qf));
650 
651   // Copy name
652   CeedCall(CeedStringAllocCopy(name, (char **)&(*qf)->gallery_name));
653   (*qf)->is_gallery = true;
654   return CEED_ERROR_SUCCESS;
655 }
656 
657 /**
658   @brief Create an identity CeedQFunction.
659 
660   Inputs are written into outputs in the order given.
661   This is useful for CeedOperators that can be represented with only the action of a CeedElemRestriction and CeedBasis, such as restriction
662 and prolongation operators for p-multigrid.
663   Backends may optimize CeedOperators with this CeedQFunction to avoid the copy of input data to output fields by using the same memory location for
664 both.
665 
666   @param[in]  ceed     Ceed object where the CeedQFunction will be created
667   @param[in]  size     Size of the QFunction fields
668   @param[in]  in_mode  CeedEvalMode for input to CeedQFunction
669   @param[in]  out_mode CeedEvalMode for output to CeedQFunction
670   @param[out] qf       Address of the variable where the newly created CeedQFunction will be stored
671 
672   @return An error code: 0 - success, otherwise - failure
673 
674   @ref User
675 **/
676 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode in_mode, CeedEvalMode out_mode, CeedQFunction *qf) {
677   CeedQFunctionContext  ctx;
678   CeedContextFieldLabel size_label;
679 
680   CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Identity", qf));
681   CeedCall(CeedQFunctionAddInput(*qf, "input", size, in_mode));
682   CeedCall(CeedQFunctionAddOutput(*qf, "output", size, out_mode));
683 
684   (*qf)->is_identity = true;
685 
686   CeedCall(CeedQFunctionGetContext(*qf, &ctx));
687   CeedCall(CeedQFunctionContextGetFieldLabel(ctx, "size", &size_label));
688   CeedCall(CeedQFunctionContextSetInt32(ctx, size_label, &size));
689   return CEED_ERROR_SUCCESS;
690 }
691 
692 /**
693   @brief Copy the pointer to a CeedQFunction.
694 
695   Both pointers should be destroyed with `CeedQFunctionDestroy()`.
696 
697   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.
698         This CeedQFunction will be destroyed if `*qf_copy` is the only reference to this CeedQFunction.
699 
700   @param[in]  qf      CeedQFunction to copy reference to
701   @param[out] qf_copy Variable to store copied reference
702 
703   @return An error code: 0 - success, otherwise - failure
704 
705   @ref User
706 **/
707 int CeedQFunctionReferenceCopy(CeedQFunction qf, CeedQFunction *qf_copy) {
708   CeedCall(CeedQFunctionReference(qf));
709   CeedCall(CeedQFunctionDestroy(qf_copy));
710   *qf_copy = qf;
711   return CEED_ERROR_SUCCESS;
712 }
713 
714 /**
715   @brief Add a CeedQFunction input
716 
717   @param[in,out] qf         CeedQFunction
718   @param[in]     field_name Name of QFunction field
719   @param[in]     size       Size of QFunction field, (num_comp * 1) for @ref CEED_EVAL_NONE,
720 (num_comp * 1) for @ref CEED_EVAL_INTERP for an H^1 space or (num_comp * dim) for an H(div) or H(curl) space,
721 (num_comp * dim) for @ref CEED_EVAL_GRAD, or (num_comp * 1) for @ref CEED_EVAL_DIV, and
722 (num_comp * curl_dim) with curl_dim = 1 if dim < 3 else dim for @ref CEED_EVAL_CURL.
723   @param[in]     eval_mode  \ref CEED_EVAL_NONE to use values directly,
724                               \ref CEED_EVAL_INTERP to use interpolated values,
725                               \ref CEED_EVAL_GRAD to use gradients,
726                               \ref CEED_EVAL_DIV to use divergence,
727                               \ref CEED_EVAL_CURL to use curl.
728 
729   @return An error code: 0 - success, otherwise - failure
730 
731   @ref User
732 **/
733 int CeedQFunctionAddInput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
734   CeedCheck(!qf->is_immutable, qf->ceed, CEED_ERROR_MAJOR, "QFunction cannot be changed after set as immutable");
735   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || size == 1, qf->ceed, CEED_ERROR_DIMENSION, "CEED_EVAL_WEIGHT should have size 1");
736   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
737     CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "QFunction field names must be unique");
738   }
739   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
740     CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "QFunction field names must be unique");
741   }
742   CeedCall(CeedQFunctionFieldSet(&qf->input_fields[qf->num_input_fields], field_name, size, eval_mode));
743   qf->num_input_fields++;
744   return CEED_ERROR_SUCCESS;
745 }
746 
747 /**
748   @brief Add a CeedQFunction output
749 
750   @param[in,out] qf         CeedQFunction
751   @param[in]     field_name Name of QFunction field
752   @param[in]     size       Size of QFunction field, (num_comp * 1) for @ref CEED_EVAL_NONE,
753 (num_comp * 1) for @ref CEED_EVAL_INTERP for an H^1 space or (num_comp * dim) for an H(div) or H(curl) space,
754 (num_comp * dim) for @ref CEED_EVAL_GRAD, or (num_comp * 1) for @ref CEED_EVAL_DIV, and
755 (num_comp * curl_dim) with curl_dim = 1 if dim < 3 else dim for @ref CEED_EVAL_CURL.
756   @param[in]     eval_mode  \ref CEED_EVAL_NONE to use values directly,
757                               \ref CEED_EVAL_INTERP to use interpolated values,
758                               \ref CEED_EVAL_GRAD to use gradients,
759                               \ref CEED_EVAL_DIV to use divergence,
760                               \ref CEED_EVAL_CURL to use curl.
761 
762   @return An error code: 0 - success, otherwise - failure
763 
764   @ref User
765 **/
766 int CeedQFunctionAddOutput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
767   CeedCheck(!qf->is_immutable, qf->ceed, CEED_ERROR_MAJOR, "QFunction cannot be changed after set as immutable");
768   CeedCheck(eval_mode != CEED_EVAL_WEIGHT, qf->ceed, CEED_ERROR_DIMENSION, "Cannot create QFunction output with CEED_EVAL_WEIGHT");
769   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
770     CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "QFunction field names must be unique");
771   }
772   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
773     CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "QFunction field names must be unique");
774   }
775   CeedCall(CeedQFunctionFieldSet(&qf->output_fields[qf->num_output_fields], field_name, size, eval_mode));
776   qf->num_output_fields++;
777   return CEED_ERROR_SUCCESS;
778 }
779 
780 /**
781   @brief Get the CeedQFunctionFields of a CeedQFunction
782 
783   Note: Calling this function asserts that setup is complete and sets the CeedQFunction as immutable.
784 
785   @param[in]  qf                CeedQFunction
786   @param[out] num_input_fields  Variable to store number of input fields
787   @param[out] input_fields      Variable to store input fields
788   @param[out] num_output_fields Variable to store number of output fields
789   @param[out] output_fields     Variable to store output fields
790 
791   @return An error code: 0 - success, otherwise - failure
792 
793   @ref Advanced
794 **/
795 int CeedQFunctionGetFields(CeedQFunction qf, CeedInt *num_input_fields, CeedQFunctionField **input_fields, CeedInt *num_output_fields,
796                            CeedQFunctionField **output_fields) {
797   qf->is_immutable = true;
798   if (num_input_fields) *num_input_fields = qf->num_input_fields;
799   if (input_fields) *input_fields = qf->input_fields;
800   if (num_output_fields) *num_output_fields = qf->num_output_fields;
801   if (output_fields) *output_fields = qf->output_fields;
802   return CEED_ERROR_SUCCESS;
803 }
804 
805 /**
806   @brief Get the name of a CeedQFunctionField
807 
808   @param[in]  qf_field   CeedQFunctionField
809   @param[out] field_name Variable to store the field name
810 
811   @return An error code: 0 - success, otherwise - failure
812 
813   @ref Advanced
814 **/
815 int CeedQFunctionFieldGetName(CeedQFunctionField qf_field, char **field_name) {
816   *field_name = (char *)qf_field->field_name;
817   return CEED_ERROR_SUCCESS;
818 }
819 
820 /**
821   @brief Get the number of components of a CeedQFunctionField
822 
823   @param[in]  qf_field CeedQFunctionField
824   @param[out] size     Variable to store the size of the field
825 
826   @return An error code: 0 - success, otherwise - failure
827 
828   @ref Advanced
829 **/
830 int CeedQFunctionFieldGetSize(CeedQFunctionField qf_field, CeedInt *size) {
831   *size = qf_field->size;
832   return CEED_ERROR_SUCCESS;
833 }
834 
835 /**
836   @brief Get the CeedEvalMode of a CeedQFunctionField
837 
838   @param[in]  qf_field  CeedQFunctionField
839   @param[out] eval_mode Variable to store the field evaluation mode
840 
841   @return An error code: 0 - success, otherwise - failure
842 
843   @ref Advanced
844 **/
845 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qf_field, CeedEvalMode *eval_mode) {
846   *eval_mode = qf_field->eval_mode;
847   return CEED_ERROR_SUCCESS;
848 }
849 
850 /**
851   @brief Set global context for a CeedQFunction
852 
853   @param[in,out] qf  CeedQFunction
854   @param[in]     ctx Context data to set
855 
856   @return An error code: 0 - success, otherwise - failure
857 
858   @ref User
859 **/
860 int CeedQFunctionSetContext(CeedQFunction qf, CeedQFunctionContext ctx) {
861   CeedCall(CeedQFunctionContextDestroy(&qf->ctx));
862   qf->ctx = ctx;
863   if (ctx) CeedCall(CeedQFunctionContextReference(ctx));
864   return CEED_ERROR_SUCCESS;
865 }
866 
867 /**
868   @brief Set writability of CeedQFunctionContext when calling the `CeedQFunctionUser`.
869 
870   The default value is `is_writable == true`.
871 
872   Setting `is_writable == true` indicates the `CeedQFunctionUser` writes into the CeedQFunctionContextData and requires memory syncronization
873 after calling `CeedQFunctionApply()`.
874 
875   Setting 'is_writable == false' asserts that `CeedQFunctionUser` does not modify the CeedQFunctionContextData.
876   Violating this assertion may lead to inconsistent data.
877 
878   Setting `is_writable == false` may offer a performance improvement on GPU backends.
879 
880   @param[in,out] qf          CeedQFunction
881   @param[in]     is_writable Writability status
882 
883   @return An error code: 0 - success, otherwise - failure
884 
885   @ref User
886 **/
887 int CeedQFunctionSetContextWritable(CeedQFunction qf, bool is_writable) {
888   qf->is_context_writable = is_writable;
889   return CEED_ERROR_SUCCESS;
890 }
891 
892 /**
893   @brief Set estimated number of FLOPs per quadrature required to apply QFunction
894 
895   @param[in]  qf    QFunction to estimate FLOPs for
896   @param[out] flops FLOPs per quadrature point estimate
897 
898   @ref Backend
899 **/
900 int CeedQFunctionSetUserFlopsEstimate(CeedQFunction qf, CeedSize flops) {
901   CeedCheck(flops >= 0, qf->ceed, CEED_ERROR_INCOMPATIBLE, "Must set non-negative FLOPs estimate");
902   qf->user_flop_estimate = flops;
903   return CEED_ERROR_SUCCESS;
904 }
905 
906 /**
907   @brief View a CeedQFunction
908 
909   @param[in] qf     CeedQFunction to view
910   @param[in] stream Stream to write; typically stdout/stderr or a file
911 
912   @return Error code: 0 - success, otherwise - failure
913 
914   @ref User
915 **/
916 int CeedQFunctionView(CeedQFunction qf, FILE *stream) {
917   char *kernel_name;
918 
919   CeedCall(CeedQFunctionGetKernelName(qf, &kernel_name));
920   fprintf(stream, "%sCeedQFunction - %s\n", qf->is_gallery ? "Gallery " : "User ", qf->is_gallery ? qf->gallery_name : kernel_name);
921 
922   fprintf(stream, "  %" CeedInt_FMT " input field%s:\n", qf->num_input_fields, qf->num_input_fields > 1 ? "s" : "");
923   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
924     CeedCall(CeedQFunctionFieldView(qf->input_fields[i], i, 1, stream));
925   }
926 
927   fprintf(stream, "  %" CeedInt_FMT " output field%s:\n", qf->num_output_fields, qf->num_output_fields > 1 ? "s" : "");
928   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
929     CeedCall(CeedQFunctionFieldView(qf->output_fields[i], i, 0, stream));
930   }
931   return CEED_ERROR_SUCCESS;
932 }
933 
934 /**
935   @brief Get the Ceed associated with a CeedQFunction
936 
937   @param[in]  qf   CeedQFunction
938   @param[out] ceed Variable to store Ceed
939 
940   @return An error code: 0 - success, otherwise - failure
941 
942   @ref Advanced
943 **/
944 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) {
945   *ceed = qf->ceed;
946   return CEED_ERROR_SUCCESS;
947 }
948 
949 /**
950   @brief Apply the action of a CeedQFunction
951 
952   Note: Calling this function asserts that setup is complete and sets the CeedQFunction as immutable.
953 
954   @param[in]  qf CeedQFunction
955   @param[in]  Q  Number of quadrature points
956   @param[in]  u  Array of input CeedVectors
957   @param[out] v  Array of output CeedVectors
958 
959   @return An error code: 0 - success, otherwise - failure
960 
961   @ref User
962 **/
963 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, CeedVector *u, CeedVector *v) {
964   CeedCheck(qf->Apply, qf->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support QFunctionApply");
965   CeedCheck(Q % qf->vec_length == 0, qf->ceed, CEED_ERROR_DIMENSION,
966             "Number of quadrature points %" CeedInt_FMT " must be a multiple of %" CeedInt_FMT, Q, qf->vec_length);
967   qf->is_immutable = true;
968   CeedCall(qf->Apply(qf, Q, u, v));
969   return CEED_ERROR_SUCCESS;
970 }
971 
972 /**
973   @brief Destroy a CeedQFunction
974 
975   @param[in,out] qf CeedQFunction to destroy
976 
977   @return An error code: 0 - success, otherwise - failure
978 
979   @ref User
980 **/
981 int CeedQFunctionDestroy(CeedQFunction *qf) {
982   if (!*qf || --(*qf)->ref_count > 0) {
983     *qf = NULL;
984     return CEED_ERROR_SUCCESS;
985   }
986   // Backend destroy
987   if ((*qf)->Destroy) {
988     CeedCall((*qf)->Destroy(*qf));
989   }
990   // Free fields
991   for (CeedInt i = 0; i < (*qf)->num_input_fields; i++) {
992     CeedCall(CeedFree(&(*(*qf)->input_fields[i]).field_name));
993     CeedCall(CeedFree(&(*qf)->input_fields[i]));
994   }
995   for (CeedInt i = 0; i < (*qf)->num_output_fields; i++) {
996     CeedCall(CeedFree(&(*(*qf)->output_fields[i]).field_name));
997     CeedCall(CeedFree(&(*qf)->output_fields[i]));
998   }
999   CeedCall(CeedFree(&(*qf)->input_fields));
1000   CeedCall(CeedFree(&(*qf)->output_fields));
1001 
1002   // User context data object
1003   CeedCall(CeedQFunctionContextDestroy(&(*qf)->ctx));
1004 
1005   CeedCall(CeedFree(&(*qf)->user_source));
1006   CeedCall(CeedFree(&(*qf)->source_path));
1007   CeedCall(CeedFree(&(*qf)->gallery_name));
1008   CeedCall(CeedFree(&(*qf)->kernel_name));
1009   CeedCall(CeedDestroy(&(*qf)->ceed));
1010   CeedCall(CeedFree(qf));
1011   return CEED_ERROR_SUCCESS;
1012 }
1013 
1014 /// @}
1015