xref: /libCEED/include/ceed-impl.h (revision fe2413ff2a5f6e0d0808f191532abb277c4b8bc7)
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 <stdbool.h>
24 
25 #define CEED_INTERN CEED_EXTERN __attribute__((visibility ("hidden")))
26 
27 #define CEED_MAX_RESOURCE_LEN 1024
28 #define CEED_ALIGN 64
29 
30 // Lookup table field for backend functions
31 typedef struct {
32   const char *fname;
33   size_t offset;
34 } foffset;
35 
36 struct Ceed_private {
37   Ceed delegate;
38   int (*Error)(Ceed, const char *, int, const char *, int, const char *, va_list);
39   int (*Destroy)(Ceed);
40   int (*VecCreate)(CeedInt, CeedVector);
41   int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode,
42                                const CeedInt *, CeedElemRestriction);
43   int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode,
44                                       const CeedInt *, CeedElemRestriction);
45   int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *,
46                              const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis);
47   int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt,
48                        const CeedScalar *,
49                        const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis);
50   int (*QFunctionCreate)(CeedQFunction);
51   int (*OperatorCreate)(CeedOperator);
52   int refcount;
53   void *data;
54   foffset foffsets[25];
55 };
56 
57 struct CeedVector_private {
58   Ceed ceed;
59   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
60   int (*SetValue)(CeedVector, CeedScalar);
61   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
62   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
63   int (*RestoreArray)(CeedVector, CeedScalar **);
64   int (*RestoreArrayRead)(CeedVector, const CeedScalar **);
65   int (*Destroy)(CeedVector);
66   int refcount;
67   CeedInt length;
68   uint64_t state;
69   uint64_t numreaders;
70   void *data;
71 };
72 
73 struct CeedElemRestriction_private {
74   Ceed ceed;
75   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedTransposeMode,
76                CeedVector, CeedVector, CeedRequest *);
77   int (*Destroy)(CeedElemRestriction);
78   int refcount;
79   CeedInt nelem;    /* number of elements */
80   CeedInt elemsize; /* number of dofs per element */
81   CeedInt ndof;     /* size of the L-vector, can be used for checking for
82                       correct vector sizes */
83   CeedInt ncomp;    /* number of components */
84   CeedInt blksize;  /* number of elements in a batch */
85   CeedInt nblk;     /* number of blocks of elements */
86   void *data;       /* place for the backend to store any data */
87 };
88 
89 struct CeedBasis_private {
90   Ceed ceed;
91   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
92                const CeedScalar *,
93                CeedScalar *);
94   int (*Destroy)(CeedBasis);
95   int refcount;
96   bool tensorbasis;      /* flag for tensor basis */
97   CeedInt dim;           /* topological dimension */
98   CeedInt ncomp;         /* number of field components (1 for scalar fields) */
99   CeedInt P1d;           /* number of nodes in one dimension */
100   CeedInt Q1d;           /* number of quadrature points in one dimension */
101   CeedInt P;             /* total number of nodes */
102   CeedInt Q;             /* total number of quadrature points */
103   CeedScalar *qref1d;    /* Array of length Q1d holding the locations of
104                             quadrature points on the 1D reference element [-1, 1] */
105   CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on
106                             the reference element */
107   CeedScalar
108   *interp1d;  /* row-major matrix of shape [Q1d, P1d] expressing the values of
109                             nodal basis functions at quadrature points */
110   CeedScalar
111   *grad1d;    /* row-major matrix of shape [Q1d, P1d] matrix expressing derivatives of
112                             nodal basis functions at quadrature points */
113   void *data;            /* place for the backend to store any data */
114 };
115 
116 struct CeedQFunctionField_private {
117   const char *fieldname;
118   CeedInt ncomp;
119   CeedEvalMode emode;
120 };
121 
122 struct CeedQFunction_private {
123   Ceed ceed;
124   int (*Apply)(CeedQFunction, CeedInt, const CeedScalar *const *,
125                CeedScalar *const *);
126   int (*Destroy)(CeedQFunction);
127   int refcount;
128   CeedInt vlength;    // Number of quadrature points must be padded to a multiple of vlength
129   CeedQFunctionField *inputfields;
130   CeedQFunctionField *outputfields;
131   CeedInt numinputfields, numoutputfields;
132   int (*function)(void*, CeedInt, const CeedScalar *const*, CeedScalar *const*);
133   const char *focca;
134   void *ctx;      /* user context for function */
135   size_t ctxsize; /* size of user context; may be used to copy to a device */
136   void *data;     /* backend data */
137 };
138 
139 struct CeedOperatorField_private {
140   CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity
141   CeedTransposeMode lmode;       /// Transpose mode for lvector ordering
142   CeedBasis basis;               /// Basis or NULL for collocated fields
143   CeedVector
144   vec;                /// State vector for passive fields, NULL for active fields
145 };
146 
147 struct CeedOperator_private {
148   Ceed ceed;
149   int refcount;
150   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
151   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
152                        CeedVector, CeedRequest *);
153   int (*Destroy)(CeedOperator);
154   CeedOperatorField *inputfields;
155   CeedOperatorField *outputfields;
156   CeedInt numelements; /// Number of elements
157   CeedInt numqpoints;  /// Number of quadrature points over all elements
158   CeedInt nfields;     /// Number of fields that have been set
159   CeedQFunction qf;
160   CeedQFunction dqf;
161   CeedQFunction dqfT;
162   bool setupdone;
163   void *data;
164 };
165 
166 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int,
167                                 const char *, va_list);
168 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int,
169                                const char *, va_list);
170 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int,
171                               const char *, va_list);
172 CEED_INTERN int CeedSetErrorHandler(Ceed ceed,
173                                     int (eh)(Ceed, const char *, int, const char *,
174                                         int, const char *, va_list));
175 
176 #endif
177