xref: /libCEED/include/ceed-impl.h (revision d3677ae8c8463ed424a8c5269c03811212785cd9)
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/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 *func_name;
84   size_t offset;
85 } FOffset;
86 
87 // Lookup table field for object delegates
88 typedef struct {
89   char *obj_name;
90   Ceed delegate;
91 } ObjDelegate;
92 
93 struct Ceed_private {
94   const char *resource;
95   Ceed delegate;
96   Ceed parent;
97   ObjDelegate *obj_delegates;
98   int obj_delegate_count;
99   Ceed op_fallback_ceed, op_fallback_parent;
100   const char *op_fallback_resource;
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 (*QFunctionContextCreate)(CeedQFunctionContext);
120   int (*OperatorCreate)(CeedOperator);
121   int (*CompositeOperatorCreate)(CeedOperator);
122   int ref_count;
123   void *data;
124   bool is_debug;
125   bool is_deterministic;
126   char err_msg[CEED_MAX_RESOURCE_LEN];
127   FOffset *f_offsets;
128 };
129 
130 struct CeedVector_private {
131   Ceed ceed;
132   int (*HasValidArray)(CeedVector, bool *);
133   int (*HasBorrowedArrayOfType)(CeedVector, CeedMemType, bool *);
134   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
135   int (*SetValue)(CeedVector, CeedScalar);
136   int (*SyncArray)(CeedVector, CeedMemType);
137   int (*TakeArray)(CeedVector, CeedMemType, CeedScalar **);
138   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
139   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
140   int (*GetArrayWrite)(CeedVector, CeedMemType, CeedScalar **);
141   int (*RestoreArray)(CeedVector);
142   int (*RestoreArrayRead)(CeedVector);
143   int (*Norm)(CeedVector, CeedNormType, CeedScalar *);
144   int (*Scale)(CeedVector, CeedScalar);
145   int (*AXPY)(CeedVector, CeedScalar, CeedVector);
146   int (*PointwiseMult)(CeedVector, CeedVector, CeedVector);
147   int (*Reciprocal)(CeedVector);
148   int (*Destroy)(CeedVector);
149   int ref_count;
150   CeedInt length;
151   uint64_t state;
152   uint64_t num_readers;
153   void *data;
154 };
155 
156 struct CeedElemRestriction_private {
157   Ceed ceed;
158   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector,
159                CeedRequest *);
160   int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector,
161                     CeedVector, CeedRequest *);
162   int (*GetOffsets)(CeedElemRestriction, CeedMemType, const CeedInt **);
163   int (*Destroy)(CeedElemRestriction);
164   int ref_count;
165   CeedInt num_elem;      /* number of elements */
166   CeedInt elem_size;     /* number of nodes per element */
167   CeedInt num_comp;      /* number of components */
168   CeedInt comp_stride;   /* Component stride for L-vector ordering */
169   CeedInt l_size;        /* size of the L-vector, can be used for checking
170                               for correct vector sizes */
171   CeedInt blk_size;      /* number of elements in a batch */
172   CeedInt num_blk;       /* number of blocks of elements */
173   CeedInt *strides;      /* strides between [nodes, components, elements] */
174   CeedInt layout[3];     /* E-vector layout [nodes, components, elements] */
175   uint64_t num_readers;  /* number of instances of offset read only access */
176   void *data;            /* place for the backend to store any data */
177 };
178 
179 struct CeedBasis_private {
180   Ceed ceed;
181   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
182                CeedVector, CeedVector);
183   int (*Destroy)(CeedBasis);
184   int ref_count;
185   bool tensor_basis;       /* flag for tensor basis */
186   CeedInt dim;             /* topological dimension */
187   CeedElemTopology topo;   /* element topology */
188   CeedInt num_comp;        /* number of field components (1 for scalar fields) */
189   CeedInt P_1d;            /* number of nodes in one dimension */
190   CeedInt Q_1d;            /* number of quadrature points in one dimension */
191   CeedInt P;               /* total number of nodes */
192   CeedInt Q;               /* total number of quadrature points */
193   CeedScalar *q_ref_1d;    /* Array of length Q1d holding the locations of
194                                 quadrature points on the 1D reference
195                                 element [-1, 1] */
196   CeedScalar
197   *q_weight_1d; /* array of length Q1d holding the quadrature weights on
198                                the reference element */
199   CeedScalar
200   *interp;    /* row-major matrix of shape [Q, P] expressing the values of
201                    nodal basis functions at quadrature points */
202   CeedScalar
203   *interp_1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of
204                    nodal basis functions at quadrature points */
205   CeedScalar
206   *grad;      /* row-major matrix of shape [dim*Q, P] matrix expressing
207                    derivatives of nodal basis functions at quadrature points */
208   CeedScalar
209   *grad_1d;   /* row-major matrix of shape [Q1d, P1d] matrix expressing
210                    derivatives of nodal basis functions at quadrature points */
211   CeedTensorContract contract; /* tensor contraction object */
212   void *data;                  /* place for the backend to store any data */
213 };
214 
215 struct CeedTensorContract_private {
216   Ceed ceed;
217   int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt,
218                const CeedScalar *restrict, CeedTransposeMode, const CeedInt,
219                const CeedScalar *restrict, CeedScalar *restrict);
220   int (*Destroy)(CeedTensorContract);
221   int ref_count;
222   void *data;
223 };
224 
225 struct CeedQFunctionField_private {
226   const char *field_name;
227   CeedInt size;
228   CeedEvalMode eval_mode;
229 };
230 
231 struct CeedQFunction_private {
232   Ceed ceed;
233   int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *);
234   int (*SetCUDAUserFunction)(CeedQFunction, void *);
235   int (*SetHIPUserFunction)(CeedQFunction, void *);
236   int (*Destroy)(CeedQFunction);
237   int ref_count;
238   CeedInt vec_length; /* Number of quadrature points must be padded to a
239                            multiple of vec_length */
240   CeedQFunctionField *input_fields;
241   CeedQFunctionField *output_fields;
242   CeedInt num_input_fields, num_output_fields;
243   CeedQFunctionUser function;
244   const char *source_path;
245   const char *kernel_name;
246   const char *gallery_name;
247   bool is_gallery;
248   bool is_identity;
249   bool is_fortran;
250   bool is_immutable;
251   CeedQFunctionContext ctx; /* user context for function */
252   void *data;          /* place for the backend to store any data */
253 };
254 
255 struct CeedQFunctionContext_private {
256   Ceed ceed;
257   int ref_count;
258   int (*HasValidData)(CeedQFunctionContext, bool *);
259   int (*HasBorrowedDataOfType)(CeedQFunctionContext, CeedMemType, bool *);
260   int (*SetData)(CeedQFunctionContext, CeedMemType, CeedCopyMode, void *);
261   int (*TakeData)(CeedQFunctionContext, CeedMemType, void *);
262   int (*GetData)(CeedQFunctionContext, CeedMemType, void *);
263   int (*RestoreData)(CeedQFunctionContext);
264   int (*Destroy)(CeedQFunctionContext);
265   CeedInt num_fields;
266   CeedInt max_fields;
267   CeedQFunctionContextFieldDescription *field_descriptions;
268   uint64_t state;
269   size_t ctx_size;
270   void *data;
271 };
272 
273 /// Struct to handle the context data to use the Fortran QFunction stub
274 /// @ingroup CeedQFunction
275 struct CeedFortranContext_private {
276   CeedQFunctionContext inner_ctx;
277   void (*f)(void *ctx, int *nq,
278             const CeedScalar *u,const CeedScalar *u1,
279             const CeedScalar *u2,const CeedScalar *u3,
280             const CeedScalar *u4,const CeedScalar *u5,
281             const CeedScalar *u6,const CeedScalar *u7,
282             const CeedScalar *u8,const CeedScalar *u9,
283             const CeedScalar *u10,const CeedScalar *u11,
284             const CeedScalar *u12,const CeedScalar *u13,
285             const CeedScalar *u14,const CeedScalar *u15,
286             CeedScalar *v,CeedScalar *v1,CeedScalar *v2,
287             CeedScalar *v3,CeedScalar *v4,CeedScalar *v5,
288             CeedScalar *v6,CeedScalar *v7,CeedScalar *v8,
289             CeedScalar *v9, CeedScalar *v10,CeedScalar *v11,
290             CeedScalar *v12,CeedScalar *v13,CeedScalar *v14,
291             CeedScalar *v15, int *err);
292 };
293 typedef struct CeedFortranContext_private *CeedFortranContext;
294 
295 struct CeedOperatorField_private {
296   CeedElemRestriction elem_restr; /* Restriction from L-vector */
297   CeedBasis basis;                /* Basis or CEED_BASIS_COLLOCATED for
298                                        collocated fields */
299   CeedVector vec;                 /* State vector for passive fields or
300                                        CEED_VECTOR_NONE for no vector */
301   const char *field_name;          /* matching QFunction field name */
302 };
303 
304 struct CeedOperator_private {
305   Ceed ceed;
306   CeedOperator op_fallback;
307   CeedQFunction qf_fallback;
308   int ref_count;
309   int (*LinearAssembleQFunction)(CeedOperator, CeedVector *,
310                                  CeedElemRestriction *, CeedRequest *);
311   int (*LinearAssembleQFunctionUpdate)(CeedOperator, CeedVector,
312                                        CeedElemRestriction, CeedRequest *);
313   int (*LinearAssembleDiagonal)(CeedOperator, CeedVector, CeedRequest *);
314   int (*LinearAssembleAddDiagonal)(CeedOperator, CeedVector, CeedRequest *);
315   int (*LinearAssemblePointBlockDiagonal)(CeedOperator, CeedVector,
316                                           CeedRequest *);
317   int (*LinearAssembleAddPointBlockDiagonal)(CeedOperator, CeedVector,
318       CeedRequest *);
319   int (*LinearAssembleSymbolic)(CeedOperator, CeedInt *, CeedInt **, CeedInt **);
320   int (*LinearAssemble)(CeedOperator, CeedVector);
321   int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *);
322   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
323   int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
324   int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
325   int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
326   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
327                        CeedVector, CeedRequest *);
328   int (*Destroy)(CeedOperator);
329   CeedOperatorField *input_fields;
330   CeedOperatorField *output_fields;
331   CeedInt num_elem;   /* Number of elements */
332   CeedInt num_qpts;   /* Number of quadrature points over all elements */
333   CeedInt num_fields; /* Number of fields that have been set */
334   CeedQFunction qf;
335   CeedQFunction dqf;
336   CeedQFunction dqfT;
337   bool is_immutable;
338   bool is_interface_setup;
339   bool is_backend_setup;
340   bool is_composite;
341   bool has_restriction;
342   bool has_qf_assembled;
343   CeedVector qf_assembled;
344   CeedElemRestriction qf_assembled_rstr;
345   CeedOperator *sub_operators;
346   CeedInt num_suboperators;
347   void *data;
348 };
349 
350 #endif
351