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