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