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