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