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