xref: /libCEED/interface/ceed-qfunction.c (revision 28f1e9f5007ebb4143767e9346c12f0b912c29da)
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)
161     // LCOV_EXCL_START
162     return CeedError(NULL, 1, "No suitable gallery QFunction");
163   // LCOV_EXCL_STOP
164 
165   // Create QFunction
166   ierr = CeedQFunctionCreateInterior(ceed, qfunctions[matchidx].vlength,
167                                      qfunctions[matchidx].f,
168                                      qfunctions[matchidx].source, qf);
169   CeedChk(ierr);
170 
171   // QFunction specific setup
172   ierr = qfunctions[matchidx].init(ceed, name, *qf); CeedChk(ierr);
173 
174   return 0;
175 }
176 
177 /**
178   @brief Create an identity CeedQFunction. Inputs are written into outputs in
179            the order given. This is useful for CeedOperators that can be
180            represented with only the action of a CeedRestriction and CeedBasis,
181            such as restriction and prolongation operators for p-multigrid.
182            Backends may optimize CeedOperators with this CeedQFunction to avoid
183            the copy of input data to output fields by using the same memory
184            location for both.
185 
186   @param ceed       A Ceed object where the CeedQFunction will be created
187   @param size       Size of the qfunction fields
188   @param[out] qf    Address of the variable where the newly created
189                       CeedQFunction will be stored
190 
191   @return An error code: 0 - success, otherwise - failure
192 
193   @ref Basic
194 **/
195 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedQFunction *qf) {
196   int ierr;
197 
198   ierr = CeedQFunctionCreateInteriorByName(ceed, "Identity", qf); CeedChk(ierr);
199 
200   (*qf)->identity = 1;
201   if (size > 1) {
202     CeedInt *ctx;
203     ierr = CeedCalloc(1, &ctx); CeedChk(ierr);
204     ctx[0] = size;
205     ierr = CeedQFunctionSetContext(*qf, ctx, sizeof(ctx)); CeedChk(ierr);
206     (*qf)->inputfields[0]->size = size;
207     (*qf)->outputfields[0]->size = size;
208   }
209 
210   return 0;
211 }
212 
213 /**
214   @brief Set a CeedQFunction field, used by CeedQFunctionAddInput/Output
215 
216   @param f          CeedQFunctionField
217   @param fieldname  Name of QFunction field
218   @param size       Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or
219                       1 for CEED_EVAL_NONE and CEED_EVAL_INTERP)
220   @param emode      \ref CEED_EVAL_NONE to use values directly,
221                       \ref CEED_EVAL_INTERP to use interpolated values,
222                       \ref CEED_EVAL_GRAD to use gradients.
223 
224   @return An error code: 0 - success, otherwise - failure
225 
226   @ref Developer
227 **/
228 static int CeedQFunctionFieldSet(CeedQFunctionField *f,const char *fieldname,
229                                  CeedInt size, CeedEvalMode emode) {
230   size_t len = strlen(fieldname);
231   char *tmp;
232   int ierr;
233   ierr = CeedCalloc(1,f); CeedChk(ierr);
234 
235   ierr = CeedCalloc(len+1, &tmp); CeedChk(ierr);
236   memcpy(tmp, fieldname, len+1);
237   (*f)->fieldname = tmp;
238   (*f)->size = size;
239   (*f)->emode = emode;
240   return 0;
241 }
242 
243 /**
244   @brief Add a CeedQFunction input
245 
246   @param qf         CeedQFunction
247   @param fieldname  Name of QFunction field
248   @param size       Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or
249                       1 for CEED_EVAL_NONE and CEED_EVAL_INTERP)
250   @param emode      \ref CEED_EVAL_NONE to use values directly,
251                       \ref CEED_EVAL_INTERP to use interpolated values,
252                       \ref CEED_EVAL_GRAD to use gradients.
253 
254   @return An error code: 0 - success, otherwise - failure
255 
256   @ref Basic
257 **/
258 int CeedQFunctionAddInput(CeedQFunction qf, const char *fieldname, CeedInt size,
259                           CeedEvalMode emode) {
260   int ierr = CeedQFunctionFieldSet(&qf->inputfields[qf->numinputfields],
261                                    fieldname, size, emode);
262   CeedChk(ierr);
263   qf->numinputfields++;
264   return 0;
265 }
266 
267 /**
268   @brief Add a CeedQFunction output
269 
270   @param qf         CeedQFunction
271   @param fieldname  Name of QFunction field
272   @param size       Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or
273                       1 for CEED_EVAL_NONE and CEED_EVAL_INTERP)
274   @param emode      \ref CEED_EVAL_NONE to use values directly,
275                       \ref CEED_EVAL_INTERP to use interpolated values,
276                       \ref CEED_EVAL_GRAD to use gradients.
277 
278   @return An error code: 0 - success, otherwise - failure
279 
280   @ref Basic
281 **/
282 int CeedQFunctionAddOutput(CeedQFunction qf, const char *fieldname,
283                            CeedInt size, CeedEvalMode emode) {
284   if (emode == CEED_EVAL_WEIGHT)
285     // LCOV_EXCL_START
286     return CeedError(qf->ceed, 1, "Cannot create QFunction output with "
287                      "CEED_EVAL_WEIGHT");
288   // LCOV_EXCL_STOP
289   int ierr = CeedQFunctionFieldSet(&qf->outputfields[qf->numoutputfields],
290                                    fieldname, size, emode);
291   CeedChk(ierr);
292   qf->numoutputfields++;
293   return 0;
294 }
295 
296 /**
297   @brief Get the Ceed associated with a CeedQFunction
298 
299   @param qf              CeedQFunction
300   @param[out] ceed       Variable to store Ceed
301 
302   @return An error code: 0 - success, otherwise - failure
303 
304   @ref Advanced
305 **/
306 
307 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) {
308   *ceed = qf->ceed;
309   return 0;
310 }
311 
312 /**
313   @brief Get the vector length of a CeedQFunction
314 
315   @param qf            CeedQFunction
316   @param[out] vlength  Variable to store vector length
317 
318   @return An error code: 0 - success, otherwise - failure
319 
320   @ref Advanced
321 **/
322 
323 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vlength) {
324   *vlength = qf->vlength;
325   return 0;
326 }
327 
328 /**
329   @brief Get the number of inputs and outputs to a CeedQFunction
330 
331   @param qf              CeedQFunction
332   @param[out] numinput   Variable to store number of input fields
333   @param[out] numoutput  Variable to store number of output fields
334 
335   @return An error code: 0 - success, otherwise - failure
336 
337   @ref Advanced
338 **/
339 
340 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *numinput,
341                             CeedInt *numoutput) {
342   if (numinput) *numinput = qf->numinputfields;
343   if (numoutput) *numoutput = qf->numoutputfields;
344   return 0;
345 }
346 
347 /**
348   @brief Get the source path string for a CeedQFunction
349 
350   @param qf              CeedQFunction
351   @param[out] source     Variable to store source path string
352 
353   @return An error code: 0 - success, otherwise - failure
354 
355   @ref Advanced
356 **/
357 
358 int CeedQFunctionGetSourcePath(CeedQFunction qf, char* *source) {
359   *source = (char *) qf->sourcepath;
360   return 0;
361 }
362 
363 /**
364   @brief Get the User Function for a CeedQFunction
365 
366   @param qf              CeedQFunction
367   @param[out] f          Variable to store user function
368 
369   @return An error code: 0 - success, otherwise - failure
370 
371   @ref Advanced
372 **/
373 
374 int CeedQFunctionGetUserFunction(CeedQFunction qf, int (**f)()) {
375   *f = (int (*)())qf->function;
376   return 0;
377 }
378 
379 /**
380   @brief Get global context size for a CeedQFunction
381 
382   @param qf              CeedQFunction
383   @param[out] ctxsize    Variable to store size of context data values
384 
385   @return An error code: 0 - success, otherwise - failure
386 
387   @ref Advanced
388 **/
389 
390 int CeedQFunctionGetContextSize(CeedQFunction qf, size_t *ctxsize) {
391   if (qf->fortranstatus) {
392     fContext *fctx = qf->ctx;
393     *ctxsize = fctx->innerctxsize;
394   } else {
395     *ctxsize = qf->ctxsize;
396   }
397   return 0;
398 }
399 
400 /**
401   @brief Get global context for a CeedQFunction
402 
403   @param qf              CeedQFunction
404   @param[out] ctx        Variable to store context data values
405 
406   @return An error code: 0 - success, otherwise - failure
407 
408   @ref Advanced
409 **/
410 
411 int CeedQFunctionGetContext(CeedQFunction qf, void* *ctx) {
412   *ctx = qf->ctx;
413   return 0;
414 }
415 
416 /**
417   @brief Determine if Fortran interface was used
418 
419   @param qf                  CeedQFunction
420   @param[out] fortranstatus  Variable to store Fortran status
421 
422   @return An error code: 0 - success, otherwise - failure
423 
424   @ref Advanced
425 **/
426 
427 int CeedQFunctionGetFortranStatus(CeedQFunction qf, bool *fortranstatus) {
428   *fortranstatus = qf->fortranstatus;
429   return 0;
430 }
431 
432 /**
433   @brief Determine if QFunction is identity
434 
435   @param qf               CeedQFunction
436   @param[out] identity    Variable to store identity status
437 
438   @return An error code: 0 - success, otherwise - failure
439 
440   @ref Advanced
441 **/
442 
443 int CeedQFunctionGetIdentityStatus(CeedQFunction qf, bool *identity) {
444   *identity = qf->identity;
445   return 0;
446 }
447 
448 /**
449   @brief Get true user context for a CeedQFunction
450 
451   @param qf              CeedQFunction
452   @param[out] ctx        Variable to store context data values
453 
454   @return An error code: 0 - success, otherwise - failure
455 
456   @ref Advanced
457 **/
458 
459 int CeedQFunctionGetInnerContext(CeedQFunction qf, void* *ctx) {
460   if (qf->fortranstatus) {
461     fContext *fctx = qf->ctx;
462     *ctx = fctx->innerctx;
463   } else {
464     *ctx = qf->ctx;
465   }
466 
467 
468   return 0;
469 }
470 
471 /**
472   @brief Get backend data of a CeedQFunction
473 
474   @param qf              CeedQFunction
475   @param[out] data       Variable to store data
476 
477   @return An error code: 0 - success, otherwise - failure
478 
479   @ref Advanced
480 **/
481 
482 int CeedQFunctionGetData(CeedQFunction qf, void* *data) {
483   *data = qf->data;
484   return 0;
485 }
486 
487 /**
488   @brief Set backend data of a CeedQFunction
489 
490   @param[out] qf         CeedQFunction
491   @param data            Data to set
492 
493   @return An error code: 0 - success, otherwise - failure
494 
495   @ref Advanced
496 **/
497 
498 int CeedQFunctionSetData(CeedQFunction qf, void* *data) {
499   qf->data = *data;
500   return 0;
501 }
502 
503 /**
504   @brief Set global context for a CeedQFunction
505 
506   @param qf       CeedQFunction
507   @param ctx      Context data to set
508   @param ctxsize  Size of context data values
509 
510   @return An error code: 0 - success, otherwise - failure
511 
512   @ref Basic
513 **/
514 int CeedQFunctionSetContext(CeedQFunction qf, void *ctx, size_t ctxsize) {
515   qf->ctx = ctx;
516   qf->ctxsize = ctxsize;
517   return 0;
518 }
519 
520 /**
521   @brief Apply the action of a CeedQFunction
522 
523   @param qf      CeedQFunction
524   @param Q       Number of quadrature points
525   @param[in] u   Array of input data arrays
526   @param[out] v  Array of output data arrays
527 
528   @return An error code: 0 - success, otherwise - failure
529 
530   @ref Advanced
531 **/
532 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q,
533                        CeedVector *u, CeedVector *v) {
534   int ierr;
535   if (!qf->Apply)
536     // LCOV_EXCL_START
537     return CeedError(qf->ceed, 1, "Backend does not support QFunctionApply");
538   // LCOV_EXCL_STOP
539   if (Q % qf->vlength)
540     // LCOV_EXCL_START
541     return CeedError(qf->ceed, 2, "Number of quadrature points %d must be a "
542                      "multiple of %d", Q, qf->vlength);
543   // LCOV_EXCL_STOP
544   ierr = qf->Apply(qf, Q, u, v); CeedChk(ierr);
545   return 0;
546 }
547 
548 /**
549   @brief Get the CeedQFunctionFields of a CeedQFunction
550 
551   @param qf                 CeedQFunction
552   @param[out] inputfields   Variable to store inputfields
553   @param[out] outputfields  Variable to store outputfields
554 
555   @return An error code: 0 - success, otherwise - failure
556 
557   @ref Advanced
558 **/
559 
560 int CeedQFunctionGetFields(CeedQFunction qf,
561                            CeedQFunctionField* *inputfields,
562                            CeedQFunctionField* *outputfields) {
563   if (inputfields)
564     *inputfields = qf->inputfields;
565   if (outputfields)
566     *outputfields = qf->outputfields;
567   return 0;
568 }
569 
570 /**
571   @brief Get the name of a CeedQFunctionField
572 
573   @param qffield         CeedQFunctionField
574   @param[out] fieldname  Variable to store the field name
575 
576   @return An error code: 0 - success, otherwise - failure
577 
578   @ref Advanced
579 **/
580 
581 int CeedQFunctionFieldGetName(CeedQFunctionField qffield,
582                               char* *fieldname) {
583   *fieldname = (char *)qffield->fieldname;
584   return 0;
585 }
586 
587 /**
588   @brief Get the number of components of a CeedQFunctionField
589 
590   @param qffield    CeedQFunctionField
591   @param[out] size  Variable to store the size of the field
592 
593   @return An error code: 0 - success, otherwise - failure
594 
595   @ref Advanced
596 **/
597 
598 int CeedQFunctionFieldGetSize(CeedQFunctionField qffield, CeedInt *size) {
599   *size = qffield->size;
600   return 0;
601 }
602 
603 /**
604   @brief Get the CeedEvalMode of a CeedQFunctionField
605 
606   @param qffield         CeedQFunctionField
607   @param[out] emode      Variable to store the field evaluation mode
608 
609   @return An error code: 0 - success, otherwise - failure
610 
611   @ref Advanced
612 **/
613 
614 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qffield,
615                                   CeedEvalMode *emode) {
616   *emode = qffield->emode;
617   return 0;
618 }
619 
620 /**
621   @brief Destroy a CeedQFunction
622 
623   @param qf CeedQFunction to destroy
624 
625   @return An error code: 0 - success, otherwise - failure
626 
627   @ref Basic
628 **/
629 int CeedQFunctionDestroy(CeedQFunction *qf) {
630   int ierr;
631 
632   if (!*qf || --(*qf)->refcount > 0)
633     return 0;
634   // Backend destroy
635   if ((*qf)->Destroy) {
636     ierr = (*qf)->Destroy(*qf); CeedChk(ierr);
637   }
638   // Free fields
639   for (int i=0; i<(*qf)->numinputfields; i++) {
640     ierr = CeedFree(&(*(*qf)->inputfields[i]).fieldname); CeedChk(ierr);
641     ierr = CeedFree(&(*qf)->inputfields[i]); CeedChk(ierr);
642   }
643   for (int i=0; i<(*qf)->numoutputfields; i++) {
644     ierr = CeedFree(&(*(*qf)->outputfields[i]).fieldname); CeedChk(ierr);
645     ierr = CeedFree(&(*qf)->outputfields[i]); CeedChk(ierr);
646   }
647   ierr = CeedFree(&(*qf)->inputfields); CeedChk(ierr);
648   ierr = CeedFree(&(*qf)->outputfields); CeedChk(ierr);
649   // Free ctx if identity
650   if ((*qf)->identity) {
651     ierr = CeedFree(&(*qf)->ctx); CeedChk(ierr);
652   }
653 
654   ierr = CeedFree(&(*qf)->sourcepath); CeedChk(ierr);
655   ierr = CeedDestroy(&(*qf)->ceed); CeedChk(ierr);
656   ierr = CeedFree(qf); CeedChk(ierr);
657   return 0;
658 }
659 
660 /// @}
661