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