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