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