xref: /libCEED/interface/ceed-qfunction.c (revision 1203703b5dc87b4acbe66c9a27384ca8ad07798d)
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     Ceed ceed;
293 
294     CeedCall(CeedQFunctionGetCeed(qf, &ceed));
295     CeedCall(CeedLoadSourceToBuffer(ceed, source_path, source_buffer));
296   }
297   return CEED_ERROR_SUCCESS;
298 }
299 
300 /**
301   @brief Get the User Function for a `CeedQFunction`
302 
303   @param[in]  qf `CeedQFunction`
304   @param[out] f  Variable to store user function
305 
306   @return An error code: 0 - success, otherwise - failure
307 
308   @ref Backend
309 **/
310 int CeedQFunctionGetUserFunction(CeedQFunction qf, CeedQFunctionUser *f) {
311   *f = qf->function;
312   return CEED_ERROR_SUCCESS;
313 }
314 
315 /**
316   @brief Get global context for a `CeedQFunction`.
317 
318   Note: For `CeedQFunction` from the Fortran interface, this function will return the Fortran context `CeedQFunctionContext`.
319 
320   @param[in]  qf  CeedQFunction
321   @param[out] ctx Variable to store CeedQFunctionContext
322 
323   @return An error code: 0 - success, otherwise - failure
324 
325   @ref Backend
326 **/
327 int CeedQFunctionGetContext(CeedQFunction qf, CeedQFunctionContext *ctx) {
328   *ctx = qf->ctx;
329   return CEED_ERROR_SUCCESS;
330 }
331 
332 /**
333   @brief Get context data of a `CeedQFunction`
334 
335   @param[in]  qf       `CeedQFunction`
336   @param[in]  mem_type Memory type on which to access the data.
337                          If the backend uses a different memory type, this will perform a copy.
338   @param[out] data     Data on memory type mem_type
339 
340   @return An error code: 0 - success, otherwise - failure
341 
342   @ref Backend
343 **/
344 int CeedQFunctionGetContextData(CeedQFunction qf, CeedMemType mem_type, void *data) {
345   bool                 is_writable;
346   CeedQFunctionContext ctx;
347 
348   CeedCall(CeedQFunctionGetContext(qf, &ctx));
349   if (ctx) {
350     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
351     if (is_writable) {
352       CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data));
353     } else {
354       CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data));
355     }
356   } else {
357     *(void **)data = NULL;
358   }
359   return CEED_ERROR_SUCCESS;
360 }
361 
362 /**
363   @brief Restore context data of a `CeedQFunction`
364 
365   @param[in]     qf   `CeedQFunction`
366   @param[in,out] data Data to restore
367 
368   @return An error code: 0 - success, otherwise - failure
369 
370   @ref Backend
371 **/
372 int CeedQFunctionRestoreContextData(CeedQFunction qf, void *data) {
373   bool                 is_writable;
374   CeedQFunctionContext ctx;
375 
376   CeedCall(CeedQFunctionGetContext(qf, &ctx));
377   if (ctx) {
378     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
379     if (is_writable) {
380       CeedCall(CeedQFunctionContextRestoreData(ctx, data));
381     } else {
382       CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data));
383     }
384   }
385   *(void **)data = NULL;
386   return CEED_ERROR_SUCCESS;
387 }
388 
389 /**
390   @brief Get true user context for a `CeedQFunction`
391 
392   Note: For all `CeedQFunction` this function will return the user `CeedQFunctionContext` and not interface context `CeedQFunctionContext`, if any such object exists.
393 
394   @param[in]  qf  `CeedQFunction`
395   @param[out] ctx Variable to store `CeedQFunctionContext`
396 
397   @return An error code: 0 - success, otherwise - failure
398   @ref Backend
399 **/
400 int CeedQFunctionGetInnerContext(CeedQFunction qf, CeedQFunctionContext *ctx) {
401   CeedQFunctionContext qf_ctx;
402 
403   CeedCall(CeedQFunctionGetContext(qf, &qf_ctx));
404   if (qf->is_fortran) {
405     CeedFortranContext fortran_ctx = NULL;
406 
407     CeedCall(CeedQFunctionContextGetData(qf_ctx, CEED_MEM_HOST, &fortran_ctx));
408     *ctx = fortran_ctx->inner_ctx;
409     CeedCall(CeedQFunctionContextRestoreData(qf_ctx, (void *)&fortran_ctx));
410   } else {
411     *ctx = qf_ctx;
412   }
413   return CEED_ERROR_SUCCESS;
414 }
415 
416 /**
417   @brief Get inner context data of a `CeedQFunction`
418 
419   @param[in]  qf       `CeedQFunction`
420   @param[in]  mem_type Memory type on which to access the data.
421                          If the backend uses a different memory type, this will perform a copy.
422   @param[out] data     Data on memory type mem_type
423 
424   @return An error code: 0 - success, otherwise - failure
425 
426   @ref Backend
427 **/
428 int CeedQFunctionGetInnerContextData(CeedQFunction qf, CeedMemType mem_type, void *data) {
429   bool                 is_writable;
430   CeedQFunctionContext ctx;
431 
432   CeedCall(CeedQFunctionGetInnerContext(qf, &ctx));
433   if (ctx) {
434     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
435     if (is_writable) {
436       CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data));
437     } else {
438       CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data));
439     }
440   } else {
441     *(void **)data = NULL;
442   }
443   return CEED_ERROR_SUCCESS;
444 }
445 
446 /**
447   @brief Restore inner context data of a `CeedQFunction`
448 
449   @param[in]     qf   `CeedQFunction`
450   @param[in,out] data Data to restore
451 
452   @return An error code: 0 - success, otherwise - failure
453 
454   @ref Backend
455 **/
456 int CeedQFunctionRestoreInnerContextData(CeedQFunction qf, void *data) {
457   bool                 is_writable;
458   CeedQFunctionContext ctx;
459 
460   CeedCall(CeedQFunctionGetInnerContext(qf, &ctx));
461   if (ctx) {
462     CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable));
463     if (is_writable) {
464       CeedCall(CeedQFunctionContextRestoreData(ctx, data));
465     } else {
466       CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data));
467     }
468   }
469   *(void **)data = NULL;
470   return CEED_ERROR_SUCCESS;
471 }
472 
473 /**
474   @brief Determine if `CeedQFunction` is identity
475 
476   @param[in]  qf          `CeedQFunction`
477   @param[out] is_identity Variable to store identity status
478 
479   @return An error code: 0 - success, otherwise - failure
480 
481   @ref Backend
482 **/
483 int CeedQFunctionIsIdentity(CeedQFunction qf, bool *is_identity) {
484   *is_identity = qf->is_identity;
485   return CEED_ERROR_SUCCESS;
486 }
487 
488 /**
489   @brief Determine if `CeedQFunctionContext` is writable
490 
491   @param[in]  qf          `CeedQFunction`
492   @param[out] is_writable Variable to store context writeable status
493 
494   @return An error code: 0 - success, otherwise - failure
495 
496   @ref Backend
497 **/
498 int CeedQFunctionIsContextWritable(CeedQFunction qf, bool *is_writable) {
499   *is_writable = qf->is_context_writable;
500   return CEED_ERROR_SUCCESS;
501 }
502 
503 /**
504   @brief Get backend data of a `CeedQFunction`
505 
506   @param[in]  qf   `CeedQFunction`
507   @param[out] data Variable to store data
508 
509   @return An error code: 0 - success, otherwise - failure
510 
511   @ref Backend
512 **/
513 int CeedQFunctionGetData(CeedQFunction qf, void *data) {
514   *(void **)data = qf->data;
515   return CEED_ERROR_SUCCESS;
516 }
517 
518 /**
519   @brief Set backend data of a `CeedQFunction`
520 
521   @param[in,out] qf   `CeedQFunction`
522   @param[in]     data Data to set
523 
524   @return An error code: 0 - success, otherwise - failure
525 
526   @ref Backend
527 **/
528 int CeedQFunctionSetData(CeedQFunction qf, void *data) {
529   qf->data = data;
530   return CEED_ERROR_SUCCESS;
531 }
532 
533 /**
534   @brief Get a boolean value indicating if the `CeedQFunction` is immutable
535 
536   @param[in]  qf           `CeedOperator`
537   @param[out] is_immutable Variable to store immutability status
538 
539   @return An error code: 0 - success, otherwise - failure
540 
541   @ref Backend
542 **/
543 int CeedQFunctionIsImmutable(CeedQFunction qf, bool *is_immutable) {
544   *is_immutable = qf->is_immutable;
545   return CEED_ERROR_SUCCESS;
546 }
547 
548 /**
549   @brief Set the immutable flag of a `CeedQFunction` to `true`
550 
551   @param[in,out] qf `CeedQFunction`
552 
553   @return An error code: 0 - success, otherwise - failure
554 
555   @ref Backend
556 **/
557 int CeedQFunctionSetImmutable(CeedQFunction qf) {
558   qf->is_immutable = true;
559   return CEED_ERROR_SUCCESS;
560 }
561 
562 /**
563   @brief Increment the reference counter for a `CeedQFunction`
564 
565   @param[in,out] qf `CeedQFunction` to increment the reference counter
566 
567   @return An error code: 0 - success, otherwise - failure
568 
569   @ref Backend
570 **/
571 int CeedQFunctionReference(CeedQFunction qf) {
572   qf->ref_count++;
573   return CEED_ERROR_SUCCESS;
574 }
575 
576 /**
577   @brief Estimate number of FLOPs per quadrature required to apply `CeedQFunction`
578 
579   @param[in]  qf    `CeedQFunction` to estimate FLOPs for
580   @param[out] flops Address of variable to hold FLOPs estimate
581 
582   @ref Backend
583 **/
584 int CeedQFunctionGetFlopsEstimate(CeedQFunction qf, CeedSize *flops) {
585   *flops = qf->user_flop_estimate;
586   return CEED_ERROR_SUCCESS;
587 }
588 
589 /// @}
590 
591 /// ----------------------------------------------------------------------------
592 /// CeedQFunction Public API
593 /// ----------------------------------------------------------------------------
594 /// @addtogroup CeedQFunctionUser
595 /// @{
596 
597 /**
598   @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms
599 
600   @param[in]  ceed       `Ceed` object used to create the `CeedQFunction`
601   @param[in]  vec_length Vector length.
602                            Caller must ensure that number of quadrature points is a multiple of `vec_length`.
603   @param[in]  f          Function pointer to evaluate action at quadrature points.
604                            See `CeedQFunctionUser`.
605   @param[in]  source     Absolute path to source of `CeedQFunctionUser`, "\abs_path\file.h:function_name".
606                            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.).
607                            The entire contents of this file and all locally included files are used during JiT compilation for GPU backends.
608                            All source files must be at the provided filepath at runtime for JiT to function.
609   @param[out] qf         Address of the variable where the newly created `CeedQFunction` will be stored
610 
611   @return An error code: 0 - success, otherwise - failure
612 
613   See \ref CeedQFunctionUser for details on the call-back function `f` arguments.
614 
615   @ref User
616 **/
617 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vec_length, CeedQFunctionUser f, const char *source, CeedQFunction *qf) {
618   char *user_source_copy;
619 
620   if (!ceed->QFunctionCreate) {
621     Ceed delegate;
622 
623     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "QFunction"));
624     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support CeedQFunctionCreateInterior");
625     CeedCall(CeedQFunctionCreateInterior(delegate, vec_length, f, source, qf));
626     return CEED_ERROR_SUCCESS;
627   }
628 
629   CeedCheck(!strlen(source) || strrchr(source, ':'), ceed, CEED_ERROR_INCOMPLETE,
630             "Provided path to source does not include function name. Provided: \"%s\"\nRequired: \"\\abs_path\\file.h:function_name\"", source);
631 
632   CeedCall(CeedCalloc(1, qf));
633   CeedCall(CeedReferenceCopy(ceed, &(*qf)->ceed));
634   (*qf)->ref_count           = 1;
635   (*qf)->vec_length          = vec_length;
636   (*qf)->is_identity         = false;
637   (*qf)->is_context_writable = true;
638   (*qf)->function            = f;
639   (*qf)->user_flop_estimate  = -1;
640   if (strlen(source)) {
641     size_t user_source_len = strlen(source);
642 
643     CeedCall(CeedCalloc(user_source_len + 1, &user_source_copy));
644     memcpy(user_source_copy, source, user_source_len);
645     (*qf)->user_source = user_source_copy;
646   }
647   CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->input_fields));
648   CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->output_fields));
649   CeedCall(ceed->QFunctionCreate(*qf));
650   return CEED_ERROR_SUCCESS;
651 }
652 
653 /**
654   @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms by name
655 
656   @param[in]  ceed `Ceed` object used to create the `CeedQFunction`
657   @param[in]  name Name of `CeedQFunction` to use from gallery
658   @param[out] qf   Address of the variable where the newly created `CeedQFunction` will be stored
659 
660   @return An error code: 0 - success, otherwise - failure
661 
662   @ref User
663 **/
664 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, CeedQFunction *qf) {
665   size_t match_len = 0, match_index = UINT_MAX;
666 
667   CeedCall(CeedQFunctionRegisterAll());
668   // Find matching backend
669   CeedCheck(name, ceed, CEED_ERROR_INCOMPLETE, "No CeedQFunction name provided");
670   for (size_t i = 0; i < num_qfunctions; i++) {
671     size_t      n;
672     const char *curr_name = gallery_qfunctions[i].name;
673     for (n = 0; curr_name[n] && curr_name[n] == name[n]; n++) {
674     }
675     if (n > match_len) {
676       match_len   = n;
677       match_index = i;
678     }
679   }
680   CeedCheck(match_len > 0, ceed, CEED_ERROR_UNSUPPORTED, "No suitable gallery CeedQFunction");
681 
682   // Create QFunction
683   CeedCall(CeedQFunctionCreateInterior(ceed, gallery_qfunctions[match_index].vec_length, gallery_qfunctions[match_index].f,
684                                        gallery_qfunctions[match_index].source, qf));
685 
686   // QFunction specific setup
687   CeedCall(gallery_qfunctions[match_index].init(ceed, name, *qf));
688 
689   // Copy name
690   CeedCall(CeedStringAllocCopy(name, (char **)&(*qf)->gallery_name));
691   (*qf)->is_gallery = true;
692   return CEED_ERROR_SUCCESS;
693 }
694 
695 /**
696   @brief Create an identity `CeedQFunction`.
697 
698   Inputs are written into outputs in the order given.
699   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.
700   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.
701 
702   @param[in]  ceed     `Ceed` object used to create the `CeedQFunction`
703   @param[in]  size     Size of the `CeedQFunction` fields
704   @param[in]  in_mode  @ref CeedEvalMode for input to `CeedQFunction`
705   @param[in]  out_mode @ref CeedEvalMode for output to `CeedQFunction`
706   @param[out] qf       Address of the variable where the newly created `CeedQFunction` will be stored
707 
708   @return An error code: 0 - success, otherwise - failure
709 
710   @ref User
711 **/
712 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode in_mode, CeedEvalMode out_mode, CeedQFunction *qf) {
713   CeedQFunctionContext  ctx;
714   CeedContextFieldLabel size_label;
715 
716   CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Identity", qf));
717   CeedCall(CeedQFunctionAddInput(*qf, "input", size, in_mode));
718   CeedCall(CeedQFunctionAddOutput(*qf, "output", size, out_mode));
719 
720   (*qf)->is_identity = true;
721 
722   CeedCall(CeedQFunctionGetContext(*qf, &ctx));
723   CeedCall(CeedQFunctionContextGetFieldLabel(ctx, "size", &size_label));
724   CeedCall(CeedQFunctionContextSetInt32(ctx, size_label, &size));
725   return CEED_ERROR_SUCCESS;
726 }
727 
728 /**
729   @brief Copy the pointer to a `CeedQFunction`.
730 
731   Both pointers should be destroyed with @ref CeedQFunctionDestroy().
732 
733   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`.
734         This `CeedQFunction` will be destroyed if `*qf_copy` is the only reference to this `CeedQFunction`.
735 
736   @param[in]  qf      `CeedQFunction` to copy reference to
737   @param[out] qf_copy Variable to store copied reference
738 
739   @return An error code: 0 - success, otherwise - failure
740 
741   @ref User
742 **/
743 int CeedQFunctionReferenceCopy(CeedQFunction qf, CeedQFunction *qf_copy) {
744   CeedCall(CeedQFunctionReference(qf));
745   CeedCall(CeedQFunctionDestroy(qf_copy));
746   *qf_copy = qf;
747   return CEED_ERROR_SUCCESS;
748 }
749 
750 /**
751   @brief Add a `CeedQFunction` input
752 
753   @param[in,out] qf         `CeedQFunction`
754   @param[in]     field_name Name of `CeedQFunction` field
755   @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.
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 CeedQFunctionAddInput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
767   bool is_immutable;
768   Ceed ceed;
769 
770   CeedCall(CeedQFunctionGetCeed(qf, &ceed));
771   CeedCall(CeedQFunctionIsImmutable(qf, &is_immutable));
772   CeedCheck(!is_immutable, ceed, CEED_ERROR_MAJOR, "QFunction cannot be changed after set as immutable");
773   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || size == 1, ceed, CEED_ERROR_DIMENSION, "CEED_EVAL_WEIGHT should have size 1");
774   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
775     CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), ceed, CEED_ERROR_MINOR, "CeedQFunction field names must be unique");
776   }
777   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
778     CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), ceed, CEED_ERROR_MINOR, "CeedQFunction field names must be unique");
779   }
780   CeedCall(CeedQFunctionFieldSet(&qf->input_fields[qf->num_input_fields], field_name, size, eval_mode));
781   qf->num_input_fields++;
782   return CEED_ERROR_SUCCESS;
783 }
784 
785 /**
786   @brief Add a `CeedQFunction` output
787 
788   @param[in,out] qf         `CeedQFunction`
789   @param[in]     field_name Name of `CeedQFunction` field
790   @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.
791   @param[in]     eval_mode  @ref CEED_EVAL_NONE to use values directly,
792                               @ref CEED_EVAL_INTERP to use interpolated values,
793                               @ref CEED_EVAL_GRAD to use gradients,
794                               @ref CEED_EVAL_DIV to use divergence,
795                               @ref CEED_EVAL_CURL to use curl.
796 
797   @return An error code: 0 - success, otherwise - failure
798 
799   @ref User
800 **/
801 int CeedQFunctionAddOutput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) {
802   bool is_immutable;
803   Ceed ceed;
804 
805   CeedCall(CeedQFunctionGetCeed(qf, &ceed));
806   CeedCall(CeedQFunctionIsImmutable(qf, &is_immutable));
807   CeedCheck(!is_immutable, ceed, CEED_ERROR_MAJOR, "CeedQFunction cannot be changed after set as immutable");
808   CeedCheck(eval_mode != CEED_EVAL_WEIGHT, ceed, CEED_ERROR_DIMENSION, "Cannot create CeedQFunction output with CEED_EVAL_WEIGHT");
809   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
810     CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), ceed, CEED_ERROR_MINOR, "CeedQFunction field names must be unique");
811   }
812   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
813     CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), ceed, CEED_ERROR_MINOR, "CeedQFunction field names must be unique");
814   }
815   CeedCall(CeedQFunctionFieldSet(&qf->output_fields[qf->num_output_fields], field_name, size, eval_mode));
816   qf->num_output_fields++;
817   return CEED_ERROR_SUCCESS;
818 }
819 
820 /**
821   @brief Get the `CeedQFunctionField` of a `CeedQFunction`
822 
823   Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable.
824 
825   @param[in]  qf                `CeedQFunction`
826   @param[out] num_input_fields  Variable to store number of input fields
827   @param[out] input_fields      Variable to store input fields
828   @param[out] num_output_fields Variable to store number of output fields
829   @param[out] output_fields     Variable to store output fields
830 
831   @return An error code: 0 - success, otherwise - failure
832 
833   @ref Advanced
834 **/
835 int CeedQFunctionGetFields(CeedQFunction qf, CeedInt *num_input_fields, CeedQFunctionField **input_fields, CeedInt *num_output_fields,
836                            CeedQFunctionField **output_fields) {
837   CeedCall(CeedQFunctionSetImmutable(qf));
838   if (num_input_fields) *num_input_fields = qf->num_input_fields;
839   if (input_fields) *input_fields = qf->input_fields;
840   if (num_output_fields) *num_output_fields = qf->num_output_fields;
841   if (output_fields) *output_fields = qf->output_fields;
842   return CEED_ERROR_SUCCESS;
843 }
844 
845 /**
846   @brief Get the name of a `CeedQFunctionField`
847 
848   @param[in]  qf_field   `CeedQFunctionField`
849   @param[out] field_name Variable to store the field name
850 
851   @return An error code: 0 - success, otherwise - failure
852 
853   @ref Advanced
854 **/
855 int CeedQFunctionFieldGetName(CeedQFunctionField qf_field, char **field_name) {
856   *field_name = (char *)qf_field->field_name;
857   return CEED_ERROR_SUCCESS;
858 }
859 
860 /**
861   @brief Get the number of components of a `CeedQFunctionField`
862 
863   @param[in]  qf_field `CeedQFunctionField`
864   @param[out] size     Variable to store the size of the field
865 
866   @return An error code: 0 - success, otherwise - failure
867 
868   @ref Advanced
869 **/
870 int CeedQFunctionFieldGetSize(CeedQFunctionField qf_field, CeedInt *size) {
871   *size = qf_field->size;
872   return CEED_ERROR_SUCCESS;
873 }
874 
875 /**
876   @brief Get the @ref CeedEvalMode of a `CeedQFunctionField`
877 
878   @param[in]  qf_field  `CeedQFunctionField`
879   @param[out] eval_mode Variable to store the field evaluation mode
880 
881   @return An error code: 0 - success, otherwise - failure
882 
883   @ref Advanced
884 **/
885 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qf_field, CeedEvalMode *eval_mode) {
886   *eval_mode = qf_field->eval_mode;
887   return CEED_ERROR_SUCCESS;
888 }
889 
890 /**
891   @brief Set global context for a `CeedQFunction`
892 
893   @param[in,out] qf  `CeedQFunction`
894   @param[in]     ctx Context data to set
895 
896   @return An error code: 0 - success, otherwise - failure
897 
898   @ref User
899 **/
900 int CeedQFunctionSetContext(CeedQFunction qf, CeedQFunctionContext ctx) {
901   CeedCall(CeedQFunctionContextDestroy(&qf->ctx));
902   qf->ctx = ctx;
903   if (ctx) CeedCall(CeedQFunctionContextReference(ctx));
904   return CEED_ERROR_SUCCESS;
905 }
906 
907 /**
908   @brief Set writability of `CeedQFunctionContext` when calling the `CeedQFunctionUser`.
909 
910   The default value is `is_writable == true`.
911 
912   Setting `is_writable == true` indicates the `CeedQFunctionUser` writes into the `CeedQFunctionContext` and requires memory synchronization after calling @ref CeedQFunctionApply().
913 
914   Setting 'is_writable == false' asserts that `CeedQFunctionUser` does not modify the `CeedQFunctionContext`.
915   Violating this assertion may lead to inconsistent data.
916 
917   Setting `is_writable == false` may offer a performance improvement on GPU backends.
918 
919   @param[in,out] qf          `CeedQFunction`
920   @param[in]     is_writable Boolean flag for writability status
921 
922   @return An error code: 0 - success, otherwise - failure
923 
924   @ref User
925 **/
926 int CeedQFunctionSetContextWritable(CeedQFunction qf, bool is_writable) {
927   qf->is_context_writable = is_writable;
928   return CEED_ERROR_SUCCESS;
929 }
930 
931 /**
932   @brief Set estimated number of FLOPs per quadrature required to apply `CeedQFunction`
933 
934   @param[in]  qf    `CeedQFunction` to estimate FLOPs for
935   @param[out] flops FLOPs per quadrature point estimate
936 
937   @ref Backend
938 **/
939 int CeedQFunctionSetUserFlopsEstimate(CeedQFunction qf, CeedSize flops) {
940   Ceed ceed;
941 
942   CeedCall(CeedQFunctionGetCeed(qf, &ceed));
943   CeedCheck(flops >= 0, ceed, CEED_ERROR_INCOMPATIBLE, "Must set non-negative FLOPs estimate");
944   qf->user_flop_estimate = flops;
945   return CEED_ERROR_SUCCESS;
946 }
947 
948 /**
949   @brief View a `CeedQFunction`
950 
951   @param[in] qf     `CeedQFunction` to view
952   @param[in] stream Stream to write; typically `stdout` or a file
953 
954   @return Error code: 0 - success, otherwise - failure
955 
956   @ref User
957 **/
958 int CeedQFunctionView(CeedQFunction qf, FILE *stream) {
959   char *kernel_name;
960 
961   CeedCall(CeedQFunctionGetKernelName(qf, &kernel_name));
962   fprintf(stream, "%sCeedQFunction - %s\n", qf->is_gallery ? "Gallery " : "User ", qf->is_gallery ? qf->gallery_name : kernel_name);
963 
964   fprintf(stream, "  %" CeedInt_FMT " input field%s:\n", qf->num_input_fields, qf->num_input_fields > 1 ? "s" : "");
965   for (CeedInt i = 0; i < qf->num_input_fields; i++) {
966     CeedCall(CeedQFunctionFieldView(qf->input_fields[i], i, 1, stream));
967   }
968 
969   fprintf(stream, "  %" CeedInt_FMT " output field%s:\n", qf->num_output_fields, qf->num_output_fields > 1 ? "s" : "");
970   for (CeedInt i = 0; i < qf->num_output_fields; i++) {
971     CeedCall(CeedQFunctionFieldView(qf->output_fields[i], i, 0, stream));
972   }
973   return CEED_ERROR_SUCCESS;
974 }
975 
976 /**
977   @brief Get the `Ceed` associated with a `CeedQFunction`
978 
979   @param[in]  qf   `CeedQFunction`
980   @param[out] ceed Variable to store`Ceed`
981 
982   @return An error code: 0 - success, otherwise - failure
983 
984   @ref Advanced
985 **/
986 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) {
987   *ceed = qf->ceed;
988   return CEED_ERROR_SUCCESS;
989 }
990 
991 /**
992   @brief Apply the action of a `CeedQFunction`
993 
994   Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable.
995 
996   @param[in]  qf `CeedQFunction`
997   @param[in]  Q  Number of quadrature points
998   @param[in]  u  Array of input `CeedVector`
999   @param[out] v  Array of output `CeedVector`
1000 
1001   @return An error code: 0 - success, otherwise - failure
1002 
1003   @ref User
1004 **/
1005 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, CeedVector *u, CeedVector *v) {
1006   CeedInt vec_length;
1007   Ceed    ceed;
1008 
1009   CeedCall(CeedQFunctionGetCeed(qf, &ceed));
1010   CeedCheck(qf->Apply, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support CeedQFunctionApply");
1011   CeedCall(CeedQFunctionGetVectorLength(qf, &vec_length));
1012   CeedCheck(Q % vec_length == 0, ceed, CEED_ERROR_DIMENSION, "Number of quadrature points %" CeedInt_FMT " must be a multiple of %" CeedInt_FMT, Q,
1013             qf->vec_length);
1014   CeedCall(CeedQFunctionSetImmutable(qf));
1015   CeedCall(qf->Apply(qf, Q, u, v));
1016   return CEED_ERROR_SUCCESS;
1017 }
1018 
1019 /**
1020   @brief Destroy a `CeedQFunction`
1021 
1022   @param[in,out] qf `CeedQFunction` to destroy
1023 
1024   @return An error code: 0 - success, otherwise - failure
1025 
1026   @ref User
1027 **/
1028 int CeedQFunctionDestroy(CeedQFunction *qf) {
1029   if (!*qf || --(*qf)->ref_count > 0) {
1030     *qf = NULL;
1031     return CEED_ERROR_SUCCESS;
1032   }
1033   // Backend destroy
1034   if ((*qf)->Destroy) {
1035     CeedCall((*qf)->Destroy(*qf));
1036   }
1037   // Free fields
1038   for (CeedInt i = 0; i < (*qf)->num_input_fields; i++) {
1039     CeedCall(CeedFree(&(*(*qf)->input_fields[i]).field_name));
1040     CeedCall(CeedFree(&(*qf)->input_fields[i]));
1041   }
1042   for (CeedInt i = 0; i < (*qf)->num_output_fields; i++) {
1043     CeedCall(CeedFree(&(*(*qf)->output_fields[i]).field_name));
1044     CeedCall(CeedFree(&(*qf)->output_fields[i]));
1045   }
1046   CeedCall(CeedFree(&(*qf)->input_fields));
1047   CeedCall(CeedFree(&(*qf)->output_fields));
1048 
1049   // User context data object
1050   CeedCall(CeedQFunctionContextDestroy(&(*qf)->ctx));
1051 
1052   CeedCall(CeedFree(&(*qf)->user_source));
1053   CeedCall(CeedFree(&(*qf)->source_path));
1054   CeedCall(CeedFree(&(*qf)->gallery_name));
1055   CeedCall(CeedFree(&(*qf)->kernel_name));
1056   CeedCall(CeedDestroy(&(*qf)->ceed));
1057   CeedCall(CeedFree(qf));
1058   return CEED_ERROR_SUCCESS;
1059 }
1060 
1061 /// @}
1062