xref: /libCEED/include/ceed-impl.h (revision 3f1466f8ca8785abaf1f35891274f21b4084d51f)
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   foffset *foffsets;
125 };
126 
127 struct CeedVector_private {
128   Ceed ceed;
129   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
130   int (*SetValue)(CeedVector, CeedScalar);
131   int (*SyncArray)(CeedVector, CeedMemType);
132   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
133   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
134   int (*RestoreArray)(CeedVector);
135   int (*RestoreArrayRead)(CeedVector);
136   int (*Norm)(CeedVector, CeedNormType, CeedScalar *);
137   int (*Destroy)(CeedVector);
138   int refcount;
139   CeedInt length;
140   uint64_t state;
141   uint64_t numreaders;
142   void *data;
143 };
144 
145 struct CeedElemRestriction_private {
146   Ceed ceed;
147   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector,
148                CeedRequest *);
149   int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector,
150                     CeedVector, CeedRequest *);
151   int (*GetOffsets)(CeedElemRestriction, CeedMemType, const CeedInt **);
152   int (*Destroy)(CeedElemRestriction);
153   int refcount;
154   CeedInt nelem;            /* number of elements */
155   CeedInt elemsize;         /* number of nodes per element */
156   CeedInt ncomp;            /* number of components */
157   CeedInt compstride;       /* Component stride for L-vector ordering */
158   CeedInt lsize;            /* size of the L-vector, can be used for checking
159                                  for correct vector sizes */
160   CeedInt blksize;          /* number of elements in a batch */
161   CeedInt nblk;             /* number of blocks of elements */
162   CeedInt *strides;         /* strides between [nodes, components, elements] */
163   CeedInt layout[3];        /* E-vector layout [nodes, components, elements] */
164   uint64_t numreaders;      /* number of instances of offset read only access */
165   void *data;               /* place for the backend to store any data */
166 };
167 
168 struct CeedBasis_private {
169   Ceed ceed;
170   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
171                CeedVector, CeedVector);
172   int (*Destroy)(CeedBasis);
173   int refcount;
174   bool tensorbasis;      /* flag for tensor basis */
175   CeedInt dim;           /* topological dimension */
176   CeedInt ncomp;         /* number of field components (1 for scalar fields) */
177   CeedInt P1d;           /* number of nodes in one dimension */
178   CeedInt Q1d;           /* number of quadrature points in one dimension */
179   CeedInt P;             /* total number of nodes */
180   CeedInt Q;             /* total number of quadrature points */
181   CeedScalar *qref1d;    /* Array of length Q1d holding the locations of
182                               quadrature points on the 1D reference
183                               element [-1, 1] */
184   CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on
185                               the reference element */
186   CeedScalar
187   *interp;    /* row-major matrix of shape [Q, P] expressing the values of
188                    nodal basis functions at quadrature points */
189   CeedScalar
190   *interp1d;  /* row-major matrix of shape [Q1d, P1d] expressing the values of
191                    nodal basis functions at quadrature points */
192   CeedScalar
193   *grad;      /* row-major matrix of shape [dim*Q, P] matrix expressing
194                    derivatives of nodal basis functions at quadrature points */
195   CeedScalar
196   *grad1d;    /* row-major matrix of shape [Q1d, P1d] matrix expressing
197                    derivatives of nodal basis functions at quadrature points */
198   CeedTensorContract contract; /* tensor contraction object */
199   void *data;                  /* place for the backend to store any data */
200 };
201 
202 struct CeedTensorContract_private {
203   Ceed ceed;
204   int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt,
205                const CeedScalar *restrict, CeedTransposeMode, const CeedInt,
206                const CeedScalar *restrict, CeedScalar *restrict);
207   int (*Destroy)(CeedTensorContract);
208   int refcount;
209   void *data;
210 };
211 
212 struct CeedQFunctionField_private {
213   const char *fieldname;
214   CeedInt size;
215   CeedEvalMode emode;
216 };
217 
218 struct CeedQFunction_private {
219   Ceed ceed;
220   int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *);
221   int (*Destroy)(CeedQFunction);
222   int refcount;
223   CeedInt vlength;    /* Number of quadrature points must be padded to a
224                            multiple of vlength */
225   CeedQFunctionField *inputfields;
226   CeedQFunctionField *outputfields;
227   CeedInt numinputfields, numoutputfields;
228   CeedQFunctionUser function;
229   const char *sourcepath;
230   const char *qfname;
231   bool fortranstatus;
232   bool identity;
233   void *ctx;      /* user context for function */
234   size_t ctxsize; /* size of user context; may be used to copy to a device */
235   void *data;     /* place for the backend to store any data */
236 };
237 
238 /// Struct to handle the context data to use the Fortran QFunction stub
239 /// @ingroup CeedQFunction
240 typedef struct {
241   CeedScalar *innerctx;
242   size_t innerctxsize;
243   void (*f)(void *ctx, int *nq,
244             const CeedScalar *u,const CeedScalar *u1,
245             const CeedScalar *u2,const CeedScalar *u3,
246             const CeedScalar *u4,const CeedScalar *u5,
247             const CeedScalar *u6,const CeedScalar *u7,
248             const CeedScalar *u8,const CeedScalar *u9,
249             const CeedScalar *u10,const CeedScalar *u11,
250             const CeedScalar *u12,const CeedScalar *u13,
251             const CeedScalar *u14,const CeedScalar *u15,
252             CeedScalar *v,CeedScalar *v1,CeedScalar *v2,
253             CeedScalar *v3,CeedScalar *v4,CeedScalar *v5,
254             CeedScalar *v6,CeedScalar *v7,CeedScalar *v8,
255             CeedScalar *v9, CeedScalar *v10,CeedScalar *v11,
256             CeedScalar *v12,CeedScalar *v13,CeedScalar *v14,
257             CeedScalar *v15, int *err);
258 } fContext;
259 
260 struct CeedOperatorField_private {
261   CeedElemRestriction Erestrict; /* Restriction from L-vector */
262   CeedBasis basis;               /* Basis or CEED_BASIS_COLLOCATED for
263                                       collocated fields */
264   CeedVector vec;                /* State vector for passive fields or
265                                       CEED_VECTOR_NONE for no vector */
266 };
267 
268 struct CeedOperator_private {
269   Ceed ceed;
270   CeedOperator opfallback;
271   CeedQFunction qffallback;
272   int refcount;
273   int (*LinearAssembleQFunction)(CeedOperator, CeedVector *,
274                                  CeedElemRestriction *, CeedRequest *);
275   int (*LinearAssembleDiagonal)(CeedOperator, CeedVector, CeedRequest *);
276   int (*LinearAssembleAddDiagonal)(CeedOperator, CeedVector, CeedRequest *);
277   int (*LinearAssemblePointBlockDiagonal)(CeedOperator, CeedVector,
278                                           CeedRequest *);
279   int (*LinearAssembleAddPointBlockDiagonal)(CeedOperator, CeedVector,
280       CeedRequest *);
281   int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *);
282   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
283   int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
284   int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
285   int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
286   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
287                        CeedVector, CeedRequest *);
288   int (*Destroy)(CeedOperator);
289   CeedOperatorField *inputfields;
290   CeedOperatorField *outputfields;
291   CeedInt numelements; /// Number of elements
292   CeedInt numqpoints;  /// Number of quadrature points over all elements
293   CeedInt nfields;     /// Number of fields that have been set
294   CeedQFunction qf;
295   CeedQFunction dqf;
296   CeedQFunction dqfT;
297   bool setupdone;
298   bool composite;
299   bool hasrestriction;
300   CeedOperator *suboperators;
301   CeedInt numsub;
302   void *data;
303 };
304 
305 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int,
306                                 const char *, va_list);
307 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int,
308                                const char *, va_list);
309 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int,
310                               const char *, va_list);
311 CEED_INTERN int CeedSetErrorHandler(Ceed ceed,
312                                     int (*eh)(Ceed, const char *, int,
313                                         const char *, int, const char *,
314                                         va_list));
315 
316 #endif
317