xref: /libCEED/include/ceed-impl.h (revision 81f92cf0318b9f3402d1efd37c448462967c3131)
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 /// @file
18 /// Private header for frontend components of libCEED
19 #ifndef _ceed_impl_h
20 #define _ceed_impl_h
21 
22 #include <ceed.h>
23 #include <ceed-backend.h>
24 #include <stdbool.h>
25 
26 /** @defgroup CeedUser Public API for Ceed
27     @ingroup Ceed
28 */
29 /** @defgroup CeedBackend Backend API for Ceed
30     @ingroup Ceed
31 */
32 /** @defgroup CeedDeveloper Internal library functions for Ceed
33     @ingroup Ceed
34 */
35 /** @defgroup CeedVectorUser Public API for CeedVector
36     @ingroup CeedVector
37 */
38 /** @defgroup CeedVectorBackend Backend API for CeedVector
39     @ingroup CeedVector
40 */
41 /** @defgroup CeedVectorDeveloper Internal library functions for CeedVector
42     @ingroup CeedVector
43 */
44 /** @defgroup CeedElemRestrictionUser Public API for CeedElemRestriction
45     @ingroup CeedElemRestriction
46 */
47 /** @defgroup CeedElemRestrictionBackend Backend API for CeedElemRestriction
48     @ingroup CeedElemRestriction
49 */
50 /** @defgroup CeedElemRestrictionDeveloper Internal library functions for CeedElemRestriction
51     @ingroup CeedElemRestriction
52 */
53 /** @defgroup CeedBasisUser Public API for CeedBasis
54     @ingroup CeedBasis
55 */
56 /** @defgroup CeedBasisBackend Backend API for CeedBasis
57     @ingroup CeedBasis
58 */
59 /** @defgroup CeedBasisDeveloper Internal library functions for CeedBasis
60     @ingroup CeedBasis
61 */
62 /** @defgroup CeedQFunctionUser Public API for CeedQFunction
63     @ingroup CeedQFunction
64 */
65 /** @defgroup CeedQFunctionBackend Backend API for CeedQFunction
66     @ingroup CeedQFunction
67 */
68 /** @defgroup CeedQFunctionDeveloper Internal library functions for CeedQFunction
69     @ingroup CeedQFunction
70 */
71 /** @defgroup CeedOperatorUser Public API for CeedOperator
72     @ingroup CeedOperator
73 */
74 /** @defgroup CeedOperatorBackend Backend API for CeedOperator
75     @ingroup CeedOperator
76 */
77 /** @defgroup CeedOperatorDeveloper Internal library functions for CeedOperator
78     @ingroup CeedOperator
79 */
80 
81 // Lookup table field for backend functions
82 typedef struct {
83   const char *fname;
84   size_t offset;
85 } foffset;
86 
87 // Lookup table field for object delegates
88 typedef struct {
89   char *objname;
90   Ceed delegate;
91 } objdelegate;
92 
93 struct Ceed_private {
94   const char *resource;
95   Ceed delegate;
96   Ceed parent;
97   objdelegate *objdelegates;
98   int objdelegatecount;
99   Ceed opfallbackceed, opfallbackparent;
100   const char *opfallbackresource;
101   int (*Error)(Ceed, const char *, int, const char *, int, const char *,
102                va_list);
103   int (*GetPreferredMemType)(CeedMemType *);
104   int (*Destroy)(Ceed);
105   int (*VectorCreate)(CeedInt, CeedVector);
106   int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode,
107                                const CeedInt *, CeedElemRestriction);
108   int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode,
109                                       const CeedInt *, CeedElemRestriction);
110   int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *,
111                              const CeedScalar *, const CeedScalar *,
112                              const CeedScalar *, CeedBasis);
113   int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt,
114                        const CeedScalar *,
115                        const CeedScalar *, const CeedScalar *,
116                        const CeedScalar *, CeedBasis);
117   int (*TensorContractCreate)(CeedBasis, CeedTensorContract);
118   int (*QFunctionCreate)(CeedQFunction);
119   int (*OperatorCreate)(CeedOperator);
120   int (*CompositeOperatorCreate)(CeedOperator);
121   int refcount;
122   bool isDeterministic;
123   void *data;
124   bool debug;
125   foffset *foffsets;
126 };
127 
128 struct CeedVector_private {
129   Ceed ceed;
130   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
131   int (*SetValue)(CeedVector, CeedScalar);
132   int (*SyncArray)(CeedVector, CeedMemType);
133   int (*TakeArray)(CeedVector, CeedMemType, CeedScalar **);
134   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
135   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
136   int (*RestoreArray)(CeedVector);
137   int (*RestoreArrayRead)(CeedVector);
138   int (*Norm)(CeedVector, CeedNormType, CeedScalar *);
139   int (*Destroy)(CeedVector);
140   int refcount;
141   CeedInt length;
142   uint64_t state;
143   uint64_t numreaders;
144   void *data;
145 };
146 
147 struct CeedElemRestriction_private {
148   Ceed ceed;
149   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector,
150                CeedRequest *);
151   int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector,
152                     CeedVector, CeedRequest *);
153   int (*GetOffsets)(CeedElemRestriction, CeedMemType, const CeedInt **);
154   int (*Destroy)(CeedElemRestriction);
155   int refcount;
156   CeedInt nelem;            /* number of elements */
157   CeedInt elemsize;         /* number of nodes per element */
158   CeedInt ncomp;            /* number of components */
159   CeedInt compstride;       /* Component stride for L-vector ordering */
160   CeedInt lsize;            /* size of the L-vector, can be used for checking
161                                  for correct vector sizes */
162   CeedInt blksize;          /* number of elements in a batch */
163   CeedInt nblk;             /* number of blocks of elements */
164   CeedInt *strides;         /* strides between [nodes, components, elements] */
165   CeedInt layout[3];        /* E-vector layout [nodes, components, elements] */
166   uint64_t numreaders;      /* number of instances of offset read only access */
167   void *data;               /* place for the backend to store any data */
168 };
169 
170 struct CeedBasis_private {
171   Ceed ceed;
172   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
173                CeedVector, CeedVector);
174   int (*Destroy)(CeedBasis);
175   int refcount;
176   bool tensorbasis;      /* flag for tensor basis */
177   CeedInt dim;           /* topological dimension */
178   CeedInt ncomp;         /* number of field components (1 for scalar fields) */
179   CeedInt P1d;           /* number of nodes in one dimension */
180   CeedInt Q1d;           /* number of quadrature points in one dimension */
181   CeedInt P;             /* total number of nodes */
182   CeedInt Q;             /* total number of quadrature points */
183   CeedScalar *qref1d;    /* Array of length Q1d holding the locations of
184                               quadrature points on the 1D reference
185                               element [-1, 1] */
186   CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on
187                               the reference element */
188   CeedScalar
189   *interp;    /* row-major matrix of shape [Q, P] expressing the values of
190                    nodal basis functions at quadrature points */
191   CeedScalar
192   *interp1d;  /* row-major matrix of shape [Q1d, P1d] expressing the values of
193                    nodal basis functions at quadrature points */
194   CeedScalar
195   *grad;      /* row-major matrix of shape [dim*Q, P] matrix expressing
196                    derivatives of nodal basis functions at quadrature points */
197   CeedScalar
198   *grad1d;    /* row-major matrix of shape [Q1d, P1d] matrix expressing
199                    derivatives of nodal basis functions at quadrature points */
200   CeedTensorContract contract; /* tensor contraction object */
201   void *data;                  /* place for the backend to store any data */
202 };
203 
204 struct CeedTensorContract_private {
205   Ceed ceed;
206   int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt,
207                const CeedScalar *restrict, CeedTransposeMode, const CeedInt,
208                const CeedScalar *restrict, CeedScalar *restrict);
209   int (*Destroy)(CeedTensorContract);
210   int refcount;
211   void *data;
212 };
213 
214 struct CeedQFunctionField_private {
215   const char *fieldname;
216   CeedInt size;
217   CeedEvalMode emode;
218 };
219 
220 struct CeedQFunction_private {
221   Ceed ceed;
222   int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *);
223   int (*Destroy)(CeedQFunction);
224   int refcount;
225   CeedInt vlength;    /* Number of quadrature points must be padded to a
226                            multiple of vlength */
227   CeedQFunctionField *inputfields;
228   CeedQFunctionField *outputfields;
229   CeedInt numinputfields, numoutputfields;
230   CeedQFunctionUser function;
231   const char *sourcepath;
232   const char *qfname;
233   bool fortranstatus;
234   bool identity;
235   void *ctx;      /* user context for function */
236   size_t ctxsize; /* size of user context; may be used to copy to a device */
237   void *data;     /* place for the backend to store any data */
238 };
239 
240 /// Struct to handle the context data to use the Fortran QFunction stub
241 /// @ingroup CeedQFunction
242 typedef struct {
243   CeedScalar *innerctx;
244   size_t innerctxsize;
245   void (*f)(void *ctx, int *nq,
246             const CeedScalar *u,const CeedScalar *u1,
247             const CeedScalar *u2,const CeedScalar *u3,
248             const CeedScalar *u4,const CeedScalar *u5,
249             const CeedScalar *u6,const CeedScalar *u7,
250             const CeedScalar *u8,const CeedScalar *u9,
251             const CeedScalar *u10,const CeedScalar *u11,
252             const CeedScalar *u12,const CeedScalar *u13,
253             const CeedScalar *u14,const CeedScalar *u15,
254             CeedScalar *v,CeedScalar *v1,CeedScalar *v2,
255             CeedScalar *v3,CeedScalar *v4,CeedScalar *v5,
256             CeedScalar *v6,CeedScalar *v7,CeedScalar *v8,
257             CeedScalar *v9, CeedScalar *v10,CeedScalar *v11,
258             CeedScalar *v12,CeedScalar *v13,CeedScalar *v14,
259             CeedScalar *v15, int *err);
260 } fContext;
261 
262 struct CeedOperatorField_private {
263   CeedElemRestriction Erestrict; /* Restriction from L-vector */
264   CeedBasis basis;               /* Basis or CEED_BASIS_COLLOCATED for
265                                       collocated fields */
266   CeedVector vec;                /* State vector for passive fields or
267                                       CEED_VECTOR_NONE for no vector */
268 };
269 
270 struct CeedOperator_private {
271   Ceed ceed;
272   CeedOperator opfallback;
273   CeedQFunction qffallback;
274   int refcount;
275   int (*LinearAssembleQFunction)(CeedOperator, CeedVector *,
276                                  CeedElemRestriction *, CeedRequest *);
277   int (*LinearAssembleDiagonal)(CeedOperator, CeedVector, CeedRequest *);
278   int (*LinearAssembleAddDiagonal)(CeedOperator, CeedVector, CeedRequest *);
279   int (*LinearAssemblePointBlockDiagonal)(CeedOperator, CeedVector,
280                                           CeedRequest *);
281   int (*LinearAssembleAddPointBlockDiagonal)(CeedOperator, CeedVector,
282       CeedRequest *);
283   int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *);
284   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
285   int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
286   int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
287   int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
288   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
289                        CeedVector, CeedRequest *);
290   int (*Destroy)(CeedOperator);
291   CeedOperatorField *inputfields;
292   CeedOperatorField *outputfields;
293   CeedInt numelements; /// Number of elements
294   CeedInt numqpoints;  /// Number of quadrature points over all elements
295   CeedInt nfields;     /// Number of fields that have been set
296   CeedQFunction qf;
297   CeedQFunction dqf;
298   CeedQFunction dqfT;
299   bool setupdone;
300   bool composite;
301   bool hasrestriction;
302   CeedOperator *suboperators;
303   CeedInt numsub;
304   void *data;
305 };
306 
307 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int,
308                                 const char *, va_list);
309 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int,
310                                const char *, va_list);
311 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int,
312                               const char *, va_list);
313 CEED_INTERN int CeedSetErrorHandler(Ceed ceed,
314                                     int (*eh)(Ceed, const char *, int,
315                                         const char *, int, const char *,
316                                         va_list));
317 
318 #endif
319