xref: /libCEED/rust/libceed-sys/c-src/include/ceed-impl.h (revision 713f43c3136a96a13935ed58fa5dd016f047cd4d) !
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 (*Destroy)(CeedVector);
81   int refcount;
82   CeedInt length;
83   uint64_t state;
84   uint64_t numreaders;
85   void *data;
86 };
87 
88 struct CeedElemRestriction_private {
89   Ceed ceed;
90   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedTransposeMode,
91                CeedVector, CeedVector, CeedRequest *);
92   int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode,
93                     CeedTransposeMode, CeedVector, CeedVector, CeedRequest *);
94   int (*Destroy)(CeedElemRestriction);
95   int refcount;
96   CeedInt nelem;    /* number of elements */
97   CeedInt elemsize; /* number of nodes per element */
98   CeedInt nnodes;   /* size of the L-vector, can be used for checking for
99                       correct vector sizes */
100   CeedInt ncomp;    /* number of components */
101   CeedInt blksize;  /* number of elements in a batch */
102   CeedInt nblk;     /* number of blocks of elements */
103   void *data;       /* place for the backend to store any data */
104 };
105 
106 struct CeedBasis_private {
107   Ceed ceed;
108   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
109                CeedVector, CeedVector);
110   int (*Destroy)(CeedBasis);
111   int refcount;
112   bool tensorbasis;      /* flag for tensor basis */
113   CeedInt dim;           /* topological dimension */
114   CeedInt ncomp;         /* number of field components (1 for scalar fields) */
115   CeedInt P1d;           /* number of nodes in one dimension */
116   CeedInt Q1d;           /* number of quadrature points in one dimension */
117   CeedInt P;             /* total number of nodes */
118   CeedInt Q;             /* total number of quadrature points */
119   CeedScalar *qref1d;    /* Array of length Q1d holding the locations of
120                             quadrature points on the 1D reference element [-1, 1] */
121   CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on
122                             the reference element */
123   CeedScalar
124   *interp;    /* row-major matrix of shape [Q, P] expressing the values of
125                             nodal basis functions at quadrature points */
126   CeedScalar
127   *interp1d;  /* row-major matrix of shape [Q1d, P1d] expressing the values of
128                             nodal basis functions at quadrature points */
129   CeedScalar
130   *grad;      /* row-major matrix of shape [dim*Q, P] matrix expressing derivatives of
131                             nodal basis functions at quadrature points */
132   CeedScalar
133   *grad1d;    /* row-major matrix of shape [Q1d, P1d] matrix expressing derivatives of
134                             nodal basis functions at quadrature points */
135   CeedTensorContract contract; /* tensor contraction object */
136   void *data;            /* place for the backend to store any data */
137 };
138 
139 struct CeedTensorContract_private {
140   Ceed ceed;
141   int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt,
142                const CeedScalar *restrict, CeedTransposeMode, const CeedInt,
143                const CeedScalar *restrict, CeedScalar *restrict);
144   int (*Destroy)(CeedTensorContract);
145   int refcount;
146   void *data;
147 };
148 
149 struct CeedQFunctionField_private {
150   const char *fieldname;
151   CeedInt size;
152   CeedEvalMode emode;
153 };
154 
155 struct CeedQFunction_private {
156   Ceed ceed;
157   int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *);
158   int (*Destroy)(CeedQFunction);
159   int refcount;
160   CeedInt vlength;    // Number of quadrature points must be padded to a multiple of vlength
161   CeedQFunctionField *inputfields;
162   CeedQFunctionField *outputfields;
163   CeedInt numinputfields, numoutputfields;
164   CeedQFunctionUser function;
165   const char *sourcepath;
166   const char *qfname;
167   bool fortranstatus;
168   bool identity;
169   void *ctx;      /* user context for function */
170   size_t ctxsize; /* size of user context; may be used to copy to a device */
171   void *data;     /* backend data */
172 };
173 
174 /// Struct to handle the context data to use the Fortran QFunction stub
175 /// @ingroup CeedQFunction
176 typedef struct {
177   CeedScalar *innerctx;
178   size_t innerctxsize;
179   void (*f)(void *ctx, int *nq,
180             const CeedScalar *u,const CeedScalar *u1,
181             const CeedScalar *u2,const CeedScalar *u3,
182             const CeedScalar *u4,const CeedScalar *u5,
183             const CeedScalar *u6,const CeedScalar *u7,
184             const CeedScalar *u8,const CeedScalar *u9,
185             const CeedScalar *u10,const CeedScalar *u11,
186             const CeedScalar *u12,const CeedScalar *u13,
187             const CeedScalar *u14,const CeedScalar *u15,
188             CeedScalar *v,CeedScalar *v1,CeedScalar *v2,
189             CeedScalar *v3,CeedScalar *v4,CeedScalar *v5,
190             CeedScalar *v6,CeedScalar *v7,CeedScalar *v8,
191             CeedScalar *v9, CeedScalar *v10,CeedScalar *v11,
192             CeedScalar *v12,CeedScalar *v13,CeedScalar *v14,
193             CeedScalar *v15, int *err);
194 } fContext;
195 
196 struct CeedOperatorField_private {
197   CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity
198   CeedTransposeMode lmode;       /// Transpose mode for lvector ordering
199   CeedBasis basis;               /// Basis or NULL for collocated fields
200   CeedVector
201   vec;                /// State vector for passive fields, NULL for active fields
202 };
203 
204 struct CeedOperator_private {
205   Ceed ceed;
206   CeedOperator opfallback;
207   CeedQFunction qffallback;
208   int refcount;
209   int (*AssembleLinearQFunction)(CeedOperator, CeedVector *,
210                                  CeedElemRestriction *, CeedRequest *);
211   int (*AssembleLinearDiagonal)(CeedOperator, CeedVector *, CeedRequest *);
212   int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *);
213   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
214   int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
215   int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
216   int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
217   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
218                        CeedVector, CeedRequest *);
219   int (*Destroy)(CeedOperator);
220   CeedOperatorField *inputfields;
221   CeedOperatorField *outputfields;
222   CeedInt numelements; /// Number of elements
223   CeedInt numqpoints;  /// Number of quadrature points over all elements
224   CeedInt nfields;     /// Number of fields that have been set
225   CeedQFunction qf;
226   CeedQFunction dqf;
227   CeedQFunction dqfT;
228   bool setupdone;
229   bool composite;
230   bool hasrestriction;
231   CeedOperator *suboperators;
232   CeedInt numsub;
233   void *data;
234 };
235 
236 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int,
237                                 const char *, va_list);
238 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int,
239                                const char *, va_list);
240 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int,
241                               const char *, va_list);
242 CEED_INTERN int CeedSetErrorHandler(Ceed ceed,
243                                     int (eh)(Ceed, const char *, int,
244                                         const char *, int, const char *,
245                                         va_list));
246 
247 #endif
248