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