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