xref: /libCEED/interface/ceed-qfunction.c (revision 07a0283791e13e6dae41f26317885b0585fe5c4e)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 #include <ceed-impl.h>
18 #include <ceed-backend.h>
19 #include <string.h>
20 #include <limits.h>
21 
22 /// @cond DOXYGEN_SKIP
23 static struct {
24   char name[CEED_MAX_RESOURCE_LEN];
25   char source[CEED_MAX_RESOURCE_LEN];
26   CeedInt vlength;
27   CeedQFunctionUser f;
28   int (*init)(Ceed ceed, const char *name, CeedQFunction qf);
29 } qfunctions[1024];
30 static size_t num_qfunctions;
31 /// @endcond
32 
33 /// @file
34 /// Implementation of public CeedQFunction interfaces
35 ///
36 /// @addtogroup CeedQFunction
37 /// @{
38 
39 /**
40   @brief Create a CeedQFunction for evaluating interior (volumetric) terms.
41 
42   @param ceed       A Ceed object where the CeedQFunction will be created
43   @param vlength    Vector length.  Caller must ensure that number of quadrature
44                     points is a multiple of vlength.
45   @param f          Function pointer to evaluate action at quadrature points.
46                     See \ref CeedQFunctionUser.
47   @param source     Absolute path to source of QFunction,
48                       "\abs_path\file.h:function_name"
49   @param[out] qf    Address of the variable where the newly created
50                       CeedQFunction will be stored
51 
52   @return An error code: 0 - success, otherwise - failure
53 
54   See \ref CeedQFunctionUser for details on the call-back function @a f's
55     arguments.
56 
57   @ref Basic
58 **/
59 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vlength,
60                                 CeedQFunctionUser f,
61                                 const char *source, CeedQFunction *qf) {
62   int ierr;
63   char *source_copy;
64 
65   if (!ceed->QFunctionCreate) {
66     Ceed delegate;
67     ierr = CeedGetObjectDelegate(ceed, &delegate, "QFunction"); CeedChk(ierr);
68 
69     if (!delegate)
70       // LCOV_EXCL_START
71       return CeedError(ceed, 1, "Backend does not support QFunctionCreate");
72     // LCOV_EXCL_STOP
73 
74     ierr = CeedQFunctionCreateInterior(delegate, vlength, f, source, qf);
75     CeedChk(ierr);
76     return 0;
77   }
78 
79   ierr = CeedCalloc(1,qf); CeedChk(ierr);
80   (*qf)->ceed = ceed;
81   ceed->refcount++;
82   (*qf)->refcount = 1;
83   (*qf)->vlength = vlength;
84   (*qf)->identity = 0;
85   (*qf)->function = f;
86   size_t slen = strlen(source) + 1;
87   ierr = CeedMalloc(slen, &source_copy); CeedChk(ierr);
88   memcpy(source_copy, source, slen);
89   (*qf)->sourcepath = source_copy;
90   ierr = CeedCalloc(16, &(*qf)->inputfields); CeedChk(ierr);
91   ierr = CeedCalloc(16, &(*qf)->outputfields); CeedChk(ierr);
92   ierr = ceed->QFunctionCreate(*qf); CeedChk(ierr);
93   return 0;
94 }
95 
96 /**
97   @brief Register a gallery QFunction
98 
99   @param name    Name for this backend to respond to
100   @param source  Absolute path to source of QFunction,
101                    "\path\CEED_DIR\gallery\folder\file.h:function_name"
102   @param vlength Vector length.  Caller must ensure that number of quadrature
103                    points is a multiple of vlength.
104   @param f       Function pointer to evaluate action at quadrature points.
105                    See \ref CeedQFunctionUser.
106   @param init    Initialization function called by CeedQFunctionInit() when the
107                    QFunction is selected.
108 
109   @return An error code: 0 - success, otherwise - failure
110 
111   @ref Advanced
112 **/
113 int CeedQFunctionRegister(const char *name, const char *source,
114                           CeedInt vlength, CeedQFunctionUser f,
115                           int (*init)(Ceed, const char *, CeedQFunction)) {
116   if (num_qfunctions >= sizeof(qfunctions) / sizeof(qfunctions[0]))
117     // LCOV_EXCL_START
118     return CeedError(NULL, 1, "Too many gallery QFunctions");
119   // LCOV_EXCL_STOP
120 
121   strncpy(qfunctions[num_qfunctions].name, name, CEED_MAX_RESOURCE_LEN);
122   qfunctions[num_qfunctions].name[CEED_MAX_RESOURCE_LEN-1] = 0;
123   strncpy(qfunctions[num_qfunctions].source, source, CEED_MAX_RESOURCE_LEN);
124   qfunctions[num_qfunctions].source[CEED_MAX_RESOURCE_LEN-1] = 0;
125   qfunctions[num_qfunctions].vlength = vlength;
126   qfunctions[num_qfunctions].f = f;
127   qfunctions[num_qfunctions].init = init;
128   num_qfunctions++;
129   return 0;
130 }
131 
132 /**
133   @brief Create a CeedQFunction for evaluating interior (volumetric) terms by name.
134 
135   @param ceed       A Ceed object where the CeedQFunction will be created
136   @param name       Name of QFunction to use from gallery
137   @param[out] qf    Address of the variable where the newly created
138                       CeedQFunction will be stored
139 
140   @return An error code: 0 - success, otherwise - failure
141 
142   @ref Basic
143 **/
144 int CeedQFunctionCreateInteriorByName(Ceed ceed,  const char *name,
145                                       CeedQFunction *qf) {
146   int ierr;
147   size_t matchlen = 0, matchidx = UINT_MAX;
148 
149   // Find matching backend
150   if (!name) return CeedError(NULL, 1, "No QFunction name provided");
151   for (size_t i=0; i<num_qfunctions; i++) {
152     size_t n;
153     const char *currname = qfunctions[i].name;
154     for (n = 0; currname[n] && currname[n] == name[n]; n++) {}
155     if (n > matchlen) {
156       matchlen = n;
157       matchidx = i;
158     }
159   }
160   if (!matchlen) return CeedError(NULL, 1, "No suitable gallery QFunction");
161 
162   // Create QFunction
163   ierr = CeedQFunctionCreateInterior(ceed, qfunctions[matchidx].vlength,
164                                      qfunctions[matchidx].f,
165                                      qfunctions[matchidx].source, qf);
166   CeedChk(ierr);
167 
168   // QFunction specific setup
169   ierr = qfunctions[matchidx].init(ceed, name, *qf); CeedChk(ierr);
170 
171   return 0;
172 }
173 
174 /**
175   @brief Create an identity CeedQFunction. Inputs are written into outputs in
176            the order given. This is useful for CeedOperators that can be
177            represented with only the action of a CeedRestriction and CeedBasis,
178            such as restriction and prolongation operators for p-multigrid.
179            Backends may optimize CeedOperators with this CeedQFunction to avoid
180            the copy of input data to output fields by using the same memory
181            location for both.
182 
183   @param ceed       A Ceed object where the CeedQFunction will be created
184   @param size       Size of the qfunction fields
185   @param[out] qf    Address of the variable where the newly created
186                       CeedQFunction will be stored
187 
188   @return An error code: 0 - success, otherwise - failure
189 
190   @ref Basic
191 **/
192 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedQFunction *qf) {
193   int ierr;
194 
195   ierr = CeedQFunctionCreateInteriorByName(ceed, "Identity", qf); CeedChk(ierr);
196 
197   (*qf)->identity = 1;
198   if (size > 1) {
199     CeedInt *ctx;
200     ierr = CeedCalloc(1, &ctx); CeedChk(ierr);
201     ctx[0] = size;
202     ierr = CeedQFunctionSetContext(*qf, ctx, sizeof(ctx)); CeedChk(ierr);
203   }
204 
205   return 0;
206 }
207 
208 /**
209   @brief Set a CeedQFunction field, used by CeedQFunctionAddInput/Output
210 
211   @param f          CeedQFunctionField
212   @param fieldname  Name of QFunction field
213   @param size       Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or
214                       1 for CEED_EVAL_NONE and CEED_EVAL_INTERP)
215   @param emode      \ref CEED_EVAL_NONE to use values directly,
216                       \ref CEED_EVAL_INTERP to use interpolated values,
217                       \ref CEED_EVAL_GRAD to use gradients.
218 
219   @return An error code: 0 - success, otherwise - failure
220 
221   @ref Developer
222 **/
223 static int CeedQFunctionFieldSet(CeedQFunctionField *f,const char *fieldname,
224                                  CeedInt size, CeedEvalMode emode) {
225   size_t len = strlen(fieldname);
226   char *tmp;
227   int ierr;
228   ierr = CeedCalloc(1,f); CeedChk(ierr);
229 
230   ierr = CeedCalloc(len+1, &tmp); CeedChk(ierr);
231   memcpy(tmp, fieldname, len+1);
232   (*f)->fieldname = tmp;
233   (*f)->size = size;
234   (*f)->emode = emode;
235   return 0;
236 }
237 
238 /**
239   @brief Add a CeedQFunction input
240 
241   @param qf         CeedQFunction
242   @param fieldname  Name of QFunction field
243   @param size       Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or
244                       1 for CEED_EVAL_NONE and CEED_EVAL_INTERP)
245   @param emode      \ref CEED_EVAL_NONE to use values directly,
246                       \ref CEED_EVAL_INTERP to use interpolated values,
247                       \ref CEED_EVAL_GRAD to use gradients.
248 
249   @return An error code: 0 - success, otherwise - failure
250 
251   @ref Basic
252 **/
253 int CeedQFunctionAddInput(CeedQFunction qf, const char *fieldname,
254                           CeedInt size, CeedEvalMode emode) {
255   int ierr = CeedQFunctionFieldSet(&qf->inputfields[qf->numinputfields],
256                                    fieldname, size, emode);
257   CeedChk(ierr);
258   qf->numinputfields++;
259   return 0;
260 }
261 
262 /**
263   @brief Add a CeedQFunction output
264 
265   @param qf         CeedQFunction
266   @param fieldname  Name of QFunction field
267   @param size       Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or
268                       1 for CEED_EVAL_NONE and CEED_EVAL_INTERP)
269   @param emode      \ref CEED_EVAL_NONE to use values directly,
270                       \ref CEED_EVAL_INTERP to use interpolated values,
271                       \ref CEED_EVAL_GRAD to use gradients.
272 
273   @return An error code: 0 - success, otherwise - failure
274 
275   @ref Basic
276 **/
277 int CeedQFunctionAddOutput(CeedQFunction qf, const char *fieldname,
278                            CeedInt size, CeedEvalMode emode) {
279   if (emode == CEED_EVAL_WEIGHT)
280     // LCOV_EXCL_START
281     return CeedError(qf->ceed, 1,
282                      "Cannot create QFunction output with CEED_EVAL_WEIGHT");
283   // LCOV_EXCL_STOP
284   int ierr = CeedQFunctionFieldSet(&qf->outputfields[qf->numoutputfields],
285                                    fieldname, size, emode);
286   CeedChk(ierr);
287   qf->numoutputfields++;
288   return 0;
289 }
290 
291 /**
292   @brief Get the Ceed associated with a CeedQFunction
293 
294   @param qf              CeedQFunction
295   @param[out] ceed       Variable to store Ceed
296 
297   @return An error code: 0 - success, otherwise - failure
298 
299   @ref Advanced
300 **/
301 
302 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) {
303   *ceed = qf->ceed;
304   return 0;
305 }
306 
307 /**
308   @brief Get the vector length of a CeedQFunction
309 
310   @param qf            CeedQFunction
311   @param[out] vlength  Variable to store vector length
312 
313   @return An error code: 0 - success, otherwise - failure
314 
315   @ref Advanced
316 **/
317 
318 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vlength) {
319   *vlength = qf->vlength;
320   return 0;
321 }
322 
323 /**
324   @brief Get the number of inputs and outputs to a CeedQFunction
325 
326   @param qf              CeedQFunction
327   @param[out] numinput   Variable to store number of input fields
328   @param[out] numoutput  Variable to store number of output fields
329 
330   @return An error code: 0 - success, otherwise - failure
331 
332   @ref Advanced
333 **/
334 
335 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *numinput,
336                             CeedInt *numoutput) {
337   if (numinput) *numinput = qf->numinputfields;
338   if (numoutput) *numoutput = qf->numoutputfields;
339   return 0;
340 }
341 
342 /**
343   @brief Get the source path string for a CeedQFunction
344 
345   @param qf              CeedQFunction
346   @param[out] source     Variable to store source path string
347 
348   @return An error code: 0 - success, otherwise - failure
349 
350   @ref Advanced
351 **/
352 
353 int CeedQFunctionGetSourcePath(CeedQFunction qf, char* *source) {
354   *source = (char *) qf->sourcepath;
355   return 0;
356 }
357 
358 /**
359   @brief Get the User Function for a CeedQFunction
360 
361   @param qf              CeedQFunction
362   @param[out] f          Variable to store user function
363 
364   @return An error code: 0 - success, otherwise - failure
365 
366   @ref Advanced
367 **/
368 
369 int CeedQFunctionGetUserFunction(CeedQFunction qf, int (**f)()) {
370   *f = (int (*)())qf->function;
371   return 0;
372 }
373 
374 /**
375   @brief Get global context size for a CeedQFunction
376 
377   @param qf              CeedQFunction
378   @param[out] ctxsize    Variable to store size of context data values
379 
380   @return An error code: 0 - success, otherwise - failure
381 
382   @ref Advanced
383 **/
384 
385 int CeedQFunctionGetContextSize(CeedQFunction qf, size_t *ctxsize) {
386   if (qf->fortranstatus) {
387     fContext *fctx = qf->ctx;
388     *ctxsize = fctx->innerctxsize;
389   } else {
390     *ctxsize = qf->ctxsize;
391   }
392   return 0;
393 }
394 
395 /**
396   @brief Get global context for a CeedQFunction
397 
398   @param qf              CeedQFunction
399   @param[out] ctx        Variable to store context data values
400 
401   @return An error code: 0 - success, otherwise - failure
402 
403   @ref Advanced
404 **/
405 
406 int CeedQFunctionGetContext(CeedQFunction qf, void* *ctx) {
407   *ctx = qf->ctx;
408   return 0;
409 }
410 
411 /**
412   @brief Determine if Fortran interface was used
413 
414   @param qf                  CeedQFunction
415   @param[out] fortranstatus  Variable to store Fortran status
416 
417   @return An error code: 0 - success, otherwise - failure
418 
419   @ref Advanced
420 **/
421 
422 int CeedQFunctionGetFortranStatus(CeedQFunction qf, bool *fortranstatus) {
423   *fortranstatus = qf->fortranstatus;
424   return 0;
425 }
426 
427 /**
428   @brief Determine if QFunction is identity
429 
430   @param qf               CeedQFunction
431   @param[out] identity    Variable to store identity status
432 
433   @return An error code: 0 - success, otherwise - failure
434 
435   @ref Advanced
436 **/
437 
438 int CeedQFunctionGetIdentityStatus(CeedQFunction qf, bool *identity) {
439   *identity = qf->identity;
440   return 0;
441 }
442 
443 /**
444   @brief Get true user context for a CeedQFunction
445 
446   @param qf              CeedQFunction
447   @param[out] ctx        Variable to store context data values
448 
449   @return An error code: 0 - success, otherwise - failure
450 
451   @ref Advanced
452 **/
453 
454 int CeedQFunctionGetInnerContext(CeedQFunction qf, void* *ctx) {
455   if (qf->fortranstatus) {
456     fContext *fctx = qf->ctx;
457     *ctx = fctx->innerctx;
458   } else {
459     *ctx = qf->ctx;
460   }
461 
462 
463   return 0;
464 }
465 
466 /**
467   @brief Get backend data of a CeedQFunction
468 
469   @param qf              CeedQFunction
470   @param[out] data       Variable to store data
471 
472   @return An error code: 0 - success, otherwise - failure
473 
474   @ref Advanced
475 **/
476 
477 int CeedQFunctionGetData(CeedQFunction qf, void* *data) {
478   *data = qf->data;
479   return 0;
480 }
481 
482 /**
483   @brief Set backend data of a CeedQFunction
484 
485   @param[out] qf         CeedQFunction
486   @param data            Data to set
487 
488   @return An error code: 0 - success, otherwise - failure
489 
490   @ref Advanced
491 **/
492 
493 int CeedQFunctionSetData(CeedQFunction qf, void* *data) {
494   qf->data = *data;
495   return 0;
496 }
497 
498 /**
499   @brief Set global context for a CeedQFunction
500 
501   @param qf       CeedQFunction
502   @param ctx      Context data to set
503   @param ctxsize  Size of context data values
504 
505   @return An error code: 0 - success, otherwise - failure
506 
507   @ref Basic
508 **/
509 int CeedQFunctionSetContext(CeedQFunction qf, void *ctx, size_t ctxsize) {
510   qf->ctx = ctx;
511   qf->ctxsize = ctxsize;
512   return 0;
513 }
514 
515 /**
516   @brief Apply the action of a CeedQFunction
517 
518   @param qf      CeedQFunction
519   @param Q       Number of quadrature points
520   @param[in] u   Array of input data arrays
521   @param[out] v  Array of output data arrays
522 
523   @return An error code: 0 - success, otherwise - failure
524 
525   @ref Advanced
526 **/
527 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q,
528                        CeedVector *u, CeedVector *v) {
529   int ierr;
530   if (!qf->Apply)
531     // LCOV_EXCL_START
532     return CeedError(qf->ceed, 1, "Backend does not support QFunctionApply");
533   // LCOV_EXCL_STOP
534   if (Q % qf->vlength)
535     // LCOV_EXCL_START
536     return CeedError(qf->ceed, 2,
537                      "Number of quadrature points %d must be a multiple of %d",
538                      Q, qf->vlength);
539   // LCOV_EXCL_STOP
540   ierr = qf->Apply(qf, Q, u, v); CeedChk(ierr);
541   return 0;
542 }
543 
544 /**
545   @brief Get the CeedQFunctionFields of a CeedQFunction
546 
547   @param qf                 CeedQFunction
548   @param[out] inputfields   Variable to store inputfields
549   @param[out] outputfields  Variable to store outputfields
550 
551   @return An error code: 0 - success, otherwise - failure
552 
553   @ref Advanced
554 **/
555 
556 int CeedQFunctionGetFields(CeedQFunction qf,
557                            CeedQFunctionField* *inputfields,
558                            CeedQFunctionField* *outputfields) {
559   if (inputfields) *inputfields = qf->inputfields;
560   if (outputfields) *outputfields = qf->outputfields;
561   return 0;
562 }
563 
564 /**
565   @brief Get the name of a CeedQFunctionField
566 
567   @param qffield         CeedQFunctionField
568   @param[out] fieldname  Variable to store the field name
569 
570   @return An error code: 0 - success, otherwise - failure
571 
572   @ref Advanced
573 **/
574 
575 int CeedQFunctionFieldGetName(CeedQFunctionField qffield,
576                               char* *fieldname) {
577   *fieldname = (char *)qffield->fieldname;
578   return 0;
579 }
580 
581 /**
582   @brief Get the number of components of a CeedQFunctionField
583 
584   @param qffield    CeedQFunctionField
585   @param[out] size  Variable to store the size of the field
586 
587   @return An error code: 0 - success, otherwise - failure
588 
589   @ref Advanced
590 **/
591 
592 int CeedQFunctionFieldGetSize(CeedQFunctionField qffield, CeedInt *size) {
593   *size = qffield->size;
594   return 0;
595 }
596 
597 /**
598   @brief Get the CeedEvalMode of a CeedQFunctionField
599 
600   @param qffield         CeedQFunctionField
601   @param[out] emode      Variable to store the field evaluation mode
602 
603   @return An error code: 0 - success, otherwise - failure
604 
605   @ref Advanced
606 **/
607 
608 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qffield,
609                                   CeedEvalMode *emode) {
610   *emode = qffield->emode;
611   return 0;
612 }
613 
614 /**
615   @brief Destroy a CeedQFunction
616 
617   @param qf CeedQFunction to destroy
618 
619   @return An error code: 0 - success, otherwise - failure
620 
621   @ref Basic
622 **/
623 int CeedQFunctionDestroy(CeedQFunction *qf) {
624   int ierr;
625 
626   if (!*qf || --(*qf)->refcount > 0) return 0;
627   // Backend destroy
628   if ((*qf)->Destroy) {
629     ierr = (*qf)->Destroy(*qf); CeedChk(ierr);
630   }
631   // Free fields
632   for (int i=0; i<(*qf)->numinputfields; i++) {
633     ierr = CeedFree(&(*(*qf)->inputfields[i]).fieldname); CeedChk(ierr);
634     ierr = CeedFree(&(*qf)->inputfields[i]); CeedChk(ierr);
635   }
636   for (int i=0; i<(*qf)->numoutputfields; i++) {
637     ierr = CeedFree(&(*(*qf)->outputfields[i]).fieldname); CeedChk(ierr);
638     ierr = CeedFree(&(*qf)->outputfields[i]); CeedChk(ierr);
639   }
640   ierr = CeedFree(&(*qf)->inputfields); CeedChk(ierr);
641   ierr = CeedFree(&(*qf)->outputfields); CeedChk(ierr);
642   // Free ctx if identity
643   if ((*qf)->identity) {
644     ierr = CeedFree(&(*qf)->ctx); CeedChk(ierr);
645   }
646 
647   ierr = CeedFree(&(*qf)->sourcepath); CeedChk(ierr);
648   ierr = CeedDestroy(&(*qf)->ceed); CeedChk(ierr);
649   ierr = CeedFree(qf); CeedChk(ierr);
650   return 0;
651 }
652 
653 /// @}
654