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