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