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