xref: /libCEED/include/ceed-impl.h (revision 6eb0d8b4aff72517bac7a1ace48e04610a8fe084)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Private header for frontend components of libCEED
10 #ifndef _ceed_impl_h
11 #define _ceed_impl_h
12 
13 #include <ceed/ceed.h>
14 #include <ceed/backend.h>
15 #include <stdbool.h>
16 
17 /** @defgroup CeedUser Public API for Ceed
18     @ingroup Ceed
19 */
20 /** @defgroup CeedBackend Backend API for Ceed
21     @ingroup Ceed
22 */
23 /** @defgroup CeedDeveloper Internal library functions for Ceed
24     @ingroup Ceed
25 */
26 /** @defgroup CeedVectorUser Public API for CeedVector
27     @ingroup CeedVector
28 */
29 /** @defgroup CeedVectorBackend Backend API for CeedVector
30     @ingroup CeedVector
31 */
32 /** @defgroup CeedVectorDeveloper Internal library functions for CeedVector
33     @ingroup CeedVector
34 */
35 /** @defgroup CeedElemRestrictionUser Public API for CeedElemRestriction
36     @ingroup CeedElemRestriction
37 */
38 /** @defgroup CeedElemRestrictionBackend Backend API for CeedElemRestriction
39     @ingroup CeedElemRestriction
40 */
41 /** @defgroup CeedElemRestrictionDeveloper Internal library functions for CeedElemRestriction
42     @ingroup CeedElemRestriction
43 */
44 /** @defgroup CeedBasisUser Public API for CeedBasis
45     @ingroup CeedBasis
46 */
47 /** @defgroup CeedBasisBackend Backend API for CeedBasis
48     @ingroup CeedBasis
49 */
50 /** @defgroup CeedBasisDeveloper Internal library functions for CeedBasis
51     @ingroup CeedBasis
52 */
53 /** @defgroup CeedQFunctionUser Public API for CeedQFunction
54     @ingroup CeedQFunction
55 */
56 /** @defgroup CeedQFunctionBackend Backend API for CeedQFunction
57     @ingroup CeedQFunction
58 */
59 /** @defgroup CeedQFunctionDeveloper Internal library functions for CeedQFunction
60     @ingroup CeedQFunction
61 */
62 /** @defgroup CeedOperatorUser Public API for CeedOperator
63     @ingroup CeedOperator
64 */
65 /** @defgroup CeedOperatorBackend Backend API for CeedOperator
66     @ingroup CeedOperator
67 */
68 /** @defgroup CeedOperatorDeveloper Internal library functions for CeedOperator
69     @ingroup CeedOperator
70 */
71 
72 // Lookup table field for backend functions
73 typedef struct {
74   const char *func_name;
75   size_t offset;
76 } FOffset;
77 
78 // Lookup table field for object delegates
79 typedef struct {
80   char *obj_name;
81   Ceed delegate;
82 } ObjDelegate;
83 
84 struct Ceed_private {
85   const char *resource;
86   Ceed delegate;
87   Ceed parent;
88   ObjDelegate *obj_delegates;
89   int obj_delegate_count;
90   Ceed op_fallback_ceed, op_fallback_parent;
91   const char *op_fallback_resource;
92   const char *jit_source_root;
93   int (*Error)(Ceed, const char *, int, const char *, int, const char *,
94                va_list *);
95   int (*GetPreferredMemType)(CeedMemType *);
96   int (*Destroy)(Ceed);
97   int (*VectorCreate)(CeedSize, CeedVector);
98   int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode,
99                                const CeedInt *, CeedElemRestriction);
100   int (*ElemRestrictionCreateOriented)(CeedMemType, CeedCopyMode,
101                                        const CeedInt *, const bool *,
102                                        CeedElemRestriction);
103   int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode,
104                                       const CeedInt *, CeedElemRestriction);
105   int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *,
106                              const CeedScalar *, const CeedScalar *,
107                              const CeedScalar *, CeedBasis);
108   int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt,
109                        const CeedScalar *,
110                        const CeedScalar *, const CeedScalar *,
111                        const CeedScalar *, CeedBasis);
112   int (*BasisCreateHdiv)(CeedElemTopology, CeedInt, CeedInt, CeedInt,
113                          const CeedScalar *,
114                          const CeedScalar *, const CeedScalar *,
115                          const CeedScalar *, CeedBasis);
116   int (*TensorContractCreate)(CeedBasis, CeedTensorContract);
117   int (*QFunctionCreate)(CeedQFunction);
118   int (*QFunctionContextCreate)(CeedQFunctionContext);
119   int (*OperatorCreate)(CeedOperator);
120   int (*CompositeOperatorCreate)(CeedOperator);
121   int ref_count;
122   void *data;
123   bool is_debug;
124   bool is_deterministic;
125   char err_msg[CEED_MAX_RESOURCE_LEN];
126   FOffset *f_offsets;
127 };
128 
129 struct CeedVector_private {
130   Ceed ceed;
131   int (*HasValidArray)(CeedVector, bool *);
132   int (*HasBorrowedArrayOfType)(CeedVector, CeedMemType, bool *);
133   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
134   int (*SetValue)(CeedVector, CeedScalar);
135   int (*SyncArray)(CeedVector, CeedMemType);
136   int (*TakeArray)(CeedVector, CeedMemType, CeedScalar **);
137   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
138   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
139   int (*GetArrayWrite)(CeedVector, CeedMemType, CeedScalar **);
140   int (*RestoreArray)(CeedVector);
141   int (*RestoreArrayRead)(CeedVector);
142   int (*Norm)(CeedVector, CeedNormType, CeedScalar *);
143   int (*Scale)(CeedVector, CeedScalar);
144   int (*AXPY)(CeedVector, CeedScalar, CeedVector);
145   int (*PointwiseMult)(CeedVector, CeedVector, CeedVector);
146   int (*Reciprocal)(CeedVector);
147   int (*Destroy)(CeedVector);
148   int ref_count;
149   CeedSize length;
150   uint64_t state;
151   uint64_t num_readers;
152   void *data;
153 };
154 
155 struct CeedElemRestriction_private {
156   Ceed ceed;
157   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector,
158                CeedRequest *);
159   int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector,
160                     CeedVector, CeedRequest *);
161   int (*GetOffsets)(CeedElemRestriction, CeedMemType, const CeedInt **);
162   int (*Destroy)(CeedElemRestriction);
163   int ref_count;
164   CeedInt num_elem;      /* number of elements */
165   CeedInt elem_size;     /* number of nodes per element */
166   CeedInt num_comp;      /* number of components */
167   CeedInt comp_stride;   /* Component stride for L-vector ordering */
168   CeedSize l_size;       /* size of the L-vector, can be used for checking
169                               for correct vector sizes */
170   CeedInt blk_size;      /* number of elements in a batch */
171   CeedInt num_blk;       /* number of blocks of elements */
172   CeedInt *strides;      /* strides between [nodes, components, elements] */
173   CeedInt layout[3];     /* E-vector layout [nodes, components, elements] */
174   uint64_t num_readers;  /* number of instances of offset read only access */
175   bool is_oriented;       /* flag for oriented restriction */
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 Q_comp;        /* number of Q-vector components (1 for H^1, dim for H(div)) */
190   CeedInt P_1d;            /* number of nodes in one dimension */
191   CeedInt Q_1d;            /* number of quadrature points in one dimension */
192   CeedInt P;               /* total number of nodes */
193   CeedInt Q;               /* total number of quadrature points */
194   CeedScalar *q_ref_1d;    /* Array of length Q1d holding the locations of
195                                 quadrature points on the 1D reference
196                                 element [-1, 1] */
197   CeedScalar
198   *q_weight_1d; /* array of length Q1d holding the quadrature weights on
199                                the reference element */
200   CeedScalar
201   *interp;    /* row-major matrix of shape [Q_comp*Q, P] expressing the values of
202                    nodal basis functions at quadrature points */
203   CeedScalar
204   *interp_1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of
205                    nodal basis functions at quadrature points */
206   CeedScalar
207   *grad;      /* row-major matrix of shape [dim*Q_comp*Q, P] matrix expressing
208                    derivatives of nodal basis functions at quadrature points */
209   CeedScalar
210   *grad_1d;   /* row-major matrix of shape [Q1d, P1d] matrix expressing
211                    derivatives of nodal basis functions at quadrature points */
212   CeedTensorContract contract; /* tensor contraction object */
213   CeedInt basis_space;  /* Initialize in basis constructor
214                         with 1,2 for H^1, H(div) FE space */
215   CeedScalar *div;  /* row-major matrix of shape [Q, P] expressing
216                         the divergence of nodal basis functions
217                         at quadrature points for H(div) discretizations */
218   void *data;                  /* place for the backend to store any data */
219 };
220 
221 struct CeedTensorContract_private {
222   Ceed ceed;
223   int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt,
224                const CeedScalar *restrict, CeedTransposeMode, const CeedInt,
225                const CeedScalar *restrict, CeedScalar *restrict);
226   int (*Destroy)(CeedTensorContract);
227   int ref_count;
228   void *data;
229 };
230 
231 struct CeedQFunctionField_private {
232   const char *field_name;
233   CeedInt size;
234   CeedEvalMode eval_mode;
235 };
236 
237 struct CeedQFunction_private {
238   Ceed ceed;
239   int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *);
240   int (*SetCUDAUserFunction)(CeedQFunction, void *);
241   int (*SetHIPUserFunction)(CeedQFunction, void *);
242   int (*Destroy)(CeedQFunction);
243   int ref_count;
244   CeedInt vec_length; /* Number of quadrature points must be padded to a
245                            multiple of vec_length */
246   CeedQFunctionField *input_fields;
247   CeedQFunctionField *output_fields;
248   CeedInt num_input_fields, num_output_fields;
249   CeedQFunctionUser function;
250   CeedInt user_flop_estimate;
251   const char *source_path;
252   const char *kernel_name;
253   const char *gallery_name;
254   bool is_gallery;
255   bool is_identity;
256   bool is_fortran;
257   bool is_immutable;
258   bool is_context_writable;
259   CeedQFunctionContext ctx; /* user context for function */
260   void *data;          /* place for the backend to store any data */
261 };
262 
263 struct CeedQFunctionContext_private {
264   Ceed ceed;
265   int ref_count;
266   int (*HasValidData)(CeedQFunctionContext, bool *);
267   int (*HasBorrowedDataOfType)(CeedQFunctionContext, CeedMemType, bool *);
268   int (*SetData)(CeedQFunctionContext, CeedMemType, CeedCopyMode, void *);
269   int (*TakeData)(CeedQFunctionContext, CeedMemType, void *);
270   int (*GetData)(CeedQFunctionContext, CeedMemType, void *);
271   int (*GetDataRead)(CeedQFunctionContext, CeedMemType, void *);
272   int (*RestoreData)(CeedQFunctionContext);
273   int (*RestoreDataRead)(CeedQFunctionContext);
274   int (*Destroy)(CeedQFunctionContext);
275   CeedInt num_fields;
276   CeedInt max_fields;
277   CeedContextFieldLabel *field_labels;
278   uint64_t state;
279   uint64_t num_readers;
280   size_t ctx_size;
281   void *data;
282 };
283 
284 /// Struct to handle the context data to use the Fortran QFunction stub
285 /// @ingroup CeedQFunction
286 struct CeedFortranContext_private {
287   CeedQFunctionContext inner_ctx;
288   void (*f)(void *ctx, int *nq,
289             const CeedScalar *u,const CeedScalar *u1,
290             const CeedScalar *u2,const CeedScalar *u3,
291             const CeedScalar *u4,const CeedScalar *u5,
292             const CeedScalar *u6,const CeedScalar *u7,
293             const CeedScalar *u8,const CeedScalar *u9,
294             const CeedScalar *u10,const CeedScalar *u11,
295             const CeedScalar *u12,const CeedScalar *u13,
296             const CeedScalar *u14,const CeedScalar *u15,
297             CeedScalar *v,CeedScalar *v1,CeedScalar *v2,
298             CeedScalar *v3,CeedScalar *v4,CeedScalar *v5,
299             CeedScalar *v6,CeedScalar *v7,CeedScalar *v8,
300             CeedScalar *v9, CeedScalar *v10,CeedScalar *v11,
301             CeedScalar *v12,CeedScalar *v13,CeedScalar *v14,
302             CeedScalar *v15, int *err);
303 };
304 typedef struct CeedFortranContext_private *CeedFortranContext;
305 
306 struct CeedContextFieldLabel_private {
307   const char *name;
308   const char *description;
309   CeedContextFieldType type;
310   size_t size;
311   size_t num_values;
312   size_t offset;
313   CeedInt num_sub_labels;
314   CeedContextFieldLabel *sub_labels;
315 };
316 
317 struct CeedOperatorField_private {
318   CeedElemRestriction elem_restr; /* Restriction from L-vector */
319   CeedBasis basis;                /* Basis or CEED_BASIS_COLLOCATED for
320                                        collocated fields */
321   CeedVector vec;                 /* State vector for passive fields or
322                                        CEED_VECTOR_NONE for no vector */
323   const char *field_name;          /* matching QFunction field name */
324 };
325 
326 struct CeedQFunctionAssemblyData_private {
327   Ceed ceed;
328   int ref_count;
329   bool is_setup;
330   bool reuse_data;
331   bool needs_data_update;
332   CeedVector vec;
333   CeedElemRestriction rstr;
334 };
335 
336 struct CeedOperator_private {
337   Ceed ceed;
338   CeedOperator op_fallback;
339   CeedQFunction qf_fallback;
340   int ref_count;
341   int (*LinearAssembleQFunction)(CeedOperator, CeedVector *,
342                                  CeedElemRestriction *, CeedRequest *);
343   int (*LinearAssembleQFunctionUpdate)(CeedOperator, CeedVector,
344                                        CeedElemRestriction, CeedRequest *);
345   int (*LinearAssembleDiagonal)(CeedOperator, CeedVector, CeedRequest *);
346   int (*LinearAssembleAddDiagonal)(CeedOperator, CeedVector, CeedRequest *);
347   int (*LinearAssemblePointBlockDiagonal)(CeedOperator, CeedVector,
348                                           CeedRequest *);
349   int (*LinearAssembleAddPointBlockDiagonal)(CeedOperator, CeedVector,
350       CeedRequest *);
351   int (*LinearAssembleSymbolic)(CeedOperator, CeedSize *, CeedInt **,
352                                 CeedInt **);
353   int (*LinearAssemble)(CeedOperator, CeedVector);
354   int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *);
355   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
356   int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
357   int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
358   int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
359   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
360                        CeedVector, CeedRequest *);
361   int (*Destroy)(CeedOperator);
362   CeedOperatorField *input_fields;
363   CeedOperatorField *output_fields;
364   CeedSize input_size, output_size;
365   CeedInt num_elem;   /* Number of elements */
366   CeedInt num_qpts;   /* Number of quadrature points over all elements */
367   CeedInt num_fields; /* Number of fields that have been set */
368   CeedQFunction qf;
369   CeedQFunction dqf;
370   CeedQFunction dqfT;
371   bool is_immutable;
372   bool is_interface_setup;
373   bool is_backend_setup;
374   bool is_composite;
375   bool has_restriction;
376   CeedQFunctionAssemblyData qf_assembled;
377   CeedOperator *sub_operators;
378   CeedInt num_suboperators;
379   void *data;
380   CeedInt num_context_labels;
381   CeedInt max_context_labels;
382   CeedContextFieldLabel *context_labels;
383 };
384 
385 #endif
386