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