xref: /libCEED/include/ceed-impl.h (revision 0036de2c0ad807bd8ab213cac28d33186520da16)
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 // Lookup table field for backend functions
27 typedef struct {
28   const char *fname;
29   size_t offset;
30 } foffset;
31 
32 // Lookup table field for object delegates
33 typedef struct {
34   char *objname;
35   Ceed delegate;
36 } objdelegate;
37 
38 struct Ceed_private {
39   const char *resource;
40   Ceed delegate;
41   Ceed parent;
42   objdelegate *objdelegates;
43   int objdelegatecount;
44   Ceed opfallbackceed, opfallbackparent;
45   const char *opfallbackresource;
46   int (*Error)(Ceed, const char *, int, const char *, int, const char *,
47                va_list);
48   int (*GetPreferredMemType)(CeedMemType *);
49   int (*Destroy)(Ceed);
50   int (*VectorCreate)(CeedInt, CeedVector);
51   int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode,
52                                const CeedInt *, CeedElemRestriction);
53   int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode,
54                                       const CeedInt *, CeedElemRestriction);
55   int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *,
56                              const CeedScalar *, const CeedScalar *,
57                              const CeedScalar *, CeedBasis);
58   int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt,
59                        const CeedScalar *,
60                        const CeedScalar *, const CeedScalar *,
61                        const CeedScalar *, CeedBasis);
62   int (*TensorContractCreate)(CeedBasis, CeedTensorContract);
63   int (*QFunctionCreate)(CeedQFunction);
64   int (*OperatorCreate)(CeedOperator);
65   int (*CompositeOperatorCreate)(CeedOperator);
66   int refcount;
67   void *data;
68   foffset *foffsets;
69 };
70 
71 struct CeedVector_private {
72   Ceed ceed;
73   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
74   int (*SetValue)(CeedVector, CeedScalar);
75   int (*SyncArray)(CeedVector, CeedMemType);
76   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
77   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
78   int (*RestoreArray)(CeedVector);
79   int (*RestoreArrayRead)(CeedVector);
80   int (*Norm)(CeedVector, CeedNormType, CeedScalar *);
81   int (*Destroy)(CeedVector);
82   int refcount;
83   CeedInt length;
84   uint64_t state;
85   uint64_t numreaders;
86   void *data;
87 };
88 
89 struct CeedElemRestriction_private {
90   Ceed ceed;
91   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector,
92                CeedRequest *);
93   int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector,
94                     CeedVector, CeedRequest *);
95   int (*Destroy)(CeedElemRestriction);
96   int refcount;
97   CeedInterlaceMode imode;  /* Interlacing mode for L-vector ordering */
98   CeedInt nelem;            /* number of elements */
99   CeedInt elemsize;         /* number of nodes per element */
100   CeedInt nnodes;           /* size of the L-vector, can be used for checking
101                                  for correct vector sizes */
102   CeedInt ncomp;            /* number of components */
103   CeedInt blksize;          /* number of elements in a batch */
104   CeedInt nblk;             /* number of blocks of elements */
105   CeedInt *strides;         /* strides between [nodes, components, elements] */
106   void *data;               /* place for the backend to store any data */
107 };
108 
109 struct CeedBasis_private {
110   Ceed ceed;
111   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
112                CeedVector, CeedVector);
113   int (*Destroy)(CeedBasis);
114   int refcount;
115   bool tensorbasis;      /* flag for tensor basis */
116   CeedInt dim;           /* topological dimension */
117   CeedInt ncomp;         /* number of field components (1 for scalar fields) */
118   CeedInt P1d;           /* number of nodes in one dimension */
119   CeedInt Q1d;           /* number of quadrature points in one dimension */
120   CeedInt P;             /* total number of nodes */
121   CeedInt Q;             /* total number of quadrature points */
122   CeedScalar *qref1d;    /* Array of length Q1d holding the locations of
123                               quadrature points on the 1D reference
124                               element [-1, 1] */
125   CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on
126                               the reference element */
127   CeedScalar
128   *interp;    /* row-major matrix of shape [Q, P] expressing the values of
129                    nodal basis functions at quadrature points */
130   CeedScalar
131   *interp1d;  /* row-major matrix of shape [Q1d, P1d] expressing the values of
132                    nodal basis functions at quadrature points */
133   CeedScalar
134   *grad;      /* row-major matrix of shape [dim*Q, P] matrix expressing
135                    derivatives of nodal basis functions at quadrature points */
136   CeedScalar
137   *grad1d;    /* row-major matrix of shape [Q1d, P1d] matrix expressing
138                    derivatives of nodal basis functions at quadrature points */
139   CeedTensorContract contract; /* tensor contraction object */
140   void *data;                  /* place for the backend to store any data */
141 };
142 
143 struct CeedTensorContract_private {
144   Ceed ceed;
145   int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt,
146                const CeedScalar *restrict, CeedTransposeMode, const CeedInt,
147                const CeedScalar *restrict, CeedScalar *restrict);
148   int (*Destroy)(CeedTensorContract);
149   int refcount;
150   void *data;
151 };
152 
153 struct CeedQFunctionField_private {
154   const char *fieldname;
155   CeedInt size;
156   CeedEvalMode emode;
157 };
158 
159 struct CeedQFunction_private {
160   Ceed ceed;
161   int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *);
162   int (*Destroy)(CeedQFunction);
163   int refcount;
164   CeedInt vlength;    /* Number of quadrature points must be padded to a
165                            multiple of vlength */
166   CeedQFunctionField *inputfields;
167   CeedQFunctionField *outputfields;
168   CeedInt numinputfields, numoutputfields;
169   CeedQFunctionUser function;
170   const char *sourcepath;
171   const char *qfname;
172   bool fortranstatus;
173   bool identity;
174   void *ctx;      /* user context for function */
175   size_t ctxsize; /* size of user context; may be used to copy to a device */
176   void *data;     /* place for the backend to store any data */
177 };
178 
179 /// Struct to handle the context data to use the Fortran QFunction stub
180 /// @ingroup CeedQFunction
181 typedef struct {
182   CeedScalar *innerctx;
183   size_t innerctxsize;
184   void (*f)(void *ctx, int *nq,
185             const CeedScalar *u,const CeedScalar *u1,
186             const CeedScalar *u2,const CeedScalar *u3,
187             const CeedScalar *u4,const CeedScalar *u5,
188             const CeedScalar *u6,const CeedScalar *u7,
189             const CeedScalar *u8,const CeedScalar *u9,
190             const CeedScalar *u10,const CeedScalar *u11,
191             const CeedScalar *u12,const CeedScalar *u13,
192             const CeedScalar *u14,const CeedScalar *u15,
193             CeedScalar *v,CeedScalar *v1,CeedScalar *v2,
194             CeedScalar *v3,CeedScalar *v4,CeedScalar *v5,
195             CeedScalar *v6,CeedScalar *v7,CeedScalar *v8,
196             CeedScalar *v9, CeedScalar *v10,CeedScalar *v11,
197             CeedScalar *v12,CeedScalar *v13,CeedScalar *v14,
198             CeedScalar *v15, int *err);
199 } fContext;
200 
201 struct CeedOperatorField_private {
202   CeedElemRestriction Erestrict; /* Restriction from L-vector */
203   CeedBasis basis;               /* Basis or CEED_BASIS_COLLOCATED for
204                                       collocated fields */
205   CeedVector vec;                /* State vector for passive fields or
206                                       CEED_VECTOR_NONE for no vector */
207 };
208 
209 struct CeedOperator_private {
210   Ceed ceed;
211   CeedOperator opfallback;
212   CeedQFunction qffallback;
213   int refcount;
214   int (*AssembleLinearQFunction)(CeedOperator, CeedVector *,
215                                  CeedElemRestriction *, CeedRequest *);
216   int (*AssembleLinearDiagonal)(CeedOperator, CeedVector *, CeedRequest *);
217   int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *);
218   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
219   int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
220   int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
221   int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
222   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
223                        CeedVector, CeedRequest *);
224   int (*Destroy)(CeedOperator);
225   CeedOperatorField *inputfields;
226   CeedOperatorField *outputfields;
227   CeedInt numelements; /// Number of elements
228   CeedInt numqpoints;  /// Number of quadrature points over all elements
229   CeedInt nfields;     /// Number of fields that have been set
230   CeedQFunction qf;
231   CeedQFunction dqf;
232   CeedQFunction dqfT;
233   bool setupdone;
234   bool composite;
235   bool hasrestriction;
236   CeedOperator *suboperators;
237   CeedInt numsub;
238   void *data;
239 };
240 
241 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int,
242                                 const char *, va_list);
243 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int,
244                                const char *, va_list);
245 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int,
246                               const char *, va_list);
247 CEED_INTERN int CeedSetErrorHandler(Ceed ceed,
248                                     int (eh)(Ceed, const char *, int,
249                                         const char *, int, const char *,
250                                         va_list));
251 
252 #endif
253