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