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