xref: /libCEED/include/ceed-impl.h (revision 9e1c8ed3e375f1dbcb31c7e047c1f7152ae5a2a7)
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 /* In the next 3 functions, p has to be the address of a pointer type, i.e. p
49    has to be a pointer to a pointer. */
50 CEED_INTERN int CeedMallocArray(size_t n, size_t unit, void *p);
51 CEED_INTERN int CeedCallocArray(size_t n, size_t unit, void *p);
52 CEED_INTERN int CeedReallocArray(size_t n, size_t unit, void *p);
53 CEED_INTERN int CeedFree(void *p);
54 
55 #define CeedChk(ierr) do { if (ierr) return ierr; } while (0)
56 /* Note that CeedMalloc and CeedCalloc will, generally, return pointers with
57    different memory alignments: CeedMalloc returns pointers aligned at
58    CEED_ALIGN bytes, while CeedCalloc uses the alignment of calloc. */
59 #define CeedMalloc(n, p) CeedMallocArray((n), sizeof(**(p)), p)
60 #define CeedCalloc(n, p) CeedCallocArray((n), sizeof(**(p)), p)
61 #define CeedRealloc(n, p) CeedReallocArray((n), sizeof(**(p)), p)
62 
63 struct CeedVector_private {
64   Ceed ceed;
65   int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *);
66   int (*SetValue)(CeedVector, CeedScalar);
67   int (*GetArray)(CeedVector, CeedMemType, CeedScalar **);
68   int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **);
69   int (*RestoreArray)(CeedVector, CeedScalar **);
70   int (*RestoreArrayRead)(CeedVector, const CeedScalar **);
71   int (*Destroy)(CeedVector);
72   int refcount;
73   CeedInt length;
74   uint64_t state;
75   void *data;
76 };
77 
78 struct CeedElemRestriction_private {
79   Ceed ceed;
80   int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedTransposeMode,
81                CeedVector, CeedVector, CeedRequest *);
82   int (*Destroy)(CeedElemRestriction);
83   int refcount;
84   CeedInt nelem;    /* number of elements */
85   CeedInt elemsize; /* number of dofs per element */
86   CeedInt ndof;     /* size of the L-vector, can be used for checking for
87                       correct vector sizes */
88   CeedInt ncomp;    /* number of components */
89   CeedInt blksize;  /* number of elements in a batch */
90   CeedInt nblk;     /* number of blocks of elements */
91   void *data;       /* place for the backend to store any data */
92 };
93 
94 struct CeedBasis_private {
95   Ceed ceed;
96   int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode,
97                const CeedScalar *,
98                CeedScalar *);
99   int (*Destroy)(CeedBasis);
100   int refcount;
101   bool tensorbasis;      /* flag for tensor basis */
102   CeedInt dim;           /* topological dimension */
103   CeedInt ncomp;         /* number of field components (1 for scalar fields) */
104   CeedInt P1d;           /* number of nodes in one dimension */
105   CeedInt Q1d;           /* number of quadrature points in one dimension */
106   CeedInt P;             /* total number of nodes */
107   CeedInt Q;             /* total number of quadrature points */
108   CeedScalar *qref1d;    /* Array of length Q1d holding the locations of
109                             quadrature points on the 1D reference element [-1, 1] */
110   CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on
111                             the reference element */
112   CeedScalar
113   *interp1d;  /* row-major matrix of shape [Q1d, P1d] expressing the values of
114                             nodal basis functions at quadrature points */
115   CeedScalar
116   *grad1d;    /* row-major matrix of shape [Q1d, P1d] matrix expressing derivatives of
117                             nodal basis functions at quadrature points */
118   void *data;            /* place for the backend to store any data */
119 };
120 
121 struct CeedQFunctionField {
122   const char *fieldname;
123   CeedInt ncomp;
124   CeedEvalMode emode;
125 };
126 
127 struct CeedQFunction_private {
128   Ceed ceed;
129   int (*Apply)(CeedQFunction, CeedInt, const CeedScalar *const *,
130                CeedScalar *const *);
131   int (*Destroy)(CeedQFunction);
132   int refcount;
133   CeedInt vlength;    // Number of quadrature points must be padded to a multiple of vlength
134   struct CeedQFunctionField inputfields[16];
135   struct CeedQFunctionField outputfields[16];
136   CeedInt numinputfields, numoutputfields;
137   int (*function)(void*, CeedInt, const CeedScalar *const*, CeedScalar *const*);
138   const char *focca;
139   void *ctx;      /* user context for function */
140   size_t ctxsize; /* size of user context; may be used to copy to a device */
141   void *data;     /* backend data */
142 };
143 
144 struct CeedOperatorField {
145   CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity
146   CeedBasis basis;               /// Basis or NULL for collocated fields
147   CeedVector
148   vec;                /// State vector for passive fields, NULL for active fields
149 };
150 
151 struct CeedOperator_private {
152   Ceed ceed;
153   int refcount;
154   int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *);
155   int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector,
156                        CeedVector, CeedRequest *);
157   int (*GetQData)(CeedOperator, CeedVector *);
158   int (*Destroy)(CeedOperator);
159   struct CeedOperatorField inputfields[16];
160   struct CeedOperatorField outputfields[16];
161   CeedInt numelements; /// Number of elements
162   CeedInt numqpoints;  /// Number of quadrature points over all elements
163   CeedInt nfields;     /// Number of fields that have been set
164   CeedQFunction qf;
165   CeedQFunction dqf;
166   CeedQFunction dqfT;
167   bool setupdone;
168   void *data;
169 };
170 
171 #endif
172