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 /// @file 18 /// Private header for frontend components of libCEED 19 #ifndef _ceed_impl_h 20 #define _ceed_impl_h 21 22 #include <ceed.h> 23 #include <stdbool.h> 24 25 #define CEED_INTERN CEED_EXTERN __attribute__((visibility ("hidden"))) 26 27 #define CEED_MAX_RESOURCE_LEN 1024 28 #define CEED_ALIGN 64 29 #define CEED_NUM_BACKEND_FUNCTIONS 25 30 31 // Lookup table field for backend functions 32 typedef struct { 33 const char *fname; 34 size_t offset; 35 } foffset; 36 37 struct Ceed_private { 38 Ceed delegate; 39 int (*Error)(Ceed, const char *, int, const char *, int, const char *, va_list); 40 int (*Destroy)(Ceed); 41 int (*VecCreate)(CeedInt, CeedVector); 42 int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode, 43 const CeedInt *, CeedElemRestriction); 44 int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode, 45 const CeedInt *, CeedElemRestriction); 46 int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *, 47 const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis); 48 int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt, 49 const CeedScalar *, 50 const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis); 51 int (*QFunctionCreate)(CeedQFunction); 52 int (*OperatorCreate)(CeedOperator); 53 int refcount; 54 void *data; 55 foffset foffsets[CEED_NUM_BACKEND_FUNCTIONS]; 56 }; 57 58 struct CeedVector_private { 59 Ceed ceed; 60 int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *); 61 int (*SetValue)(CeedVector, 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 uint64_t state; 70 uint64_t numreaders; 71 void *data; 72 }; 73 74 struct CeedElemRestriction_private { 75 Ceed ceed; 76 int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedTransposeMode, 77 CeedVector, CeedVector, CeedRequest *); 78 int (*Destroy)(CeedElemRestriction); 79 int refcount; 80 CeedInt nelem; /* number of elements */ 81 CeedInt elemsize; /* number of dofs per element */ 82 CeedInt ndof; /* size of the L-vector, can be used for checking for 83 correct vector sizes */ 84 CeedInt ncomp; /* number of components */ 85 CeedInt blksize; /* number of elements in a batch */ 86 CeedInt nblk; /* number of blocks of elements */ 87 void *data; /* place for the backend to store any data */ 88 }; 89 90 struct CeedBasis_private { 91 Ceed ceed; 92 int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, 93 CeedVector, CeedVector); 94 int (*Destroy)(CeedBasis); 95 int refcount; 96 bool tensorbasis; /* flag for tensor basis */ 97 CeedInt dim; /* topological dimension */ 98 CeedInt ncomp; /* number of field components (1 for scalar fields) */ 99 CeedInt P1d; /* number of nodes in one dimension */ 100 CeedInt Q1d; /* number of quadrature points in one dimension */ 101 CeedInt P; /* total number of nodes */ 102 CeedInt Q; /* total number of quadrature points */ 103 CeedScalar *qref1d; /* Array of length Q1d holding the locations of 104 quadrature points on the 1D reference element [-1, 1] */ 105 CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on 106 the reference element */ 107 CeedScalar 108 *interp1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of 109 nodal basis functions at quadrature points */ 110 CeedScalar 111 *grad1d; /* row-major matrix of shape [Q1d, P1d] matrix expressing derivatives of 112 nodal basis functions at quadrature points */ 113 void *data; /* place for the backend to store any data */ 114 }; 115 116 struct CeedQFunctionField_private { 117 const char *fieldname; 118 CeedInt ncomp; 119 CeedEvalMode emode; 120 }; 121 122 struct CeedQFunction_private { 123 Ceed ceed; 124 int (*Apply)(CeedQFunction, CeedInt, CeedVector *, 125 CeedVector *); 126 int (*Destroy)(CeedQFunction); 127 int refcount; 128 CeedInt vlength; // Number of quadrature points must be padded to a multiple of vlength 129 CeedQFunctionField *inputfields; 130 CeedQFunctionField *outputfields; 131 CeedInt numinputfields, numoutputfields; 132 int (*function)(void*, CeedInt, const CeedScalar *const*, CeedScalar *const*); 133 const char *focca; 134 bool fortranstatus; 135 void *ctx; /* user context for function */ 136 size_t ctxsize; /* size of user context; may be used to copy to a device */ 137 void *data; /* backend data */ 138 }; 139 140 typedef struct { 141 void *innerctx; 142 size_t innerctxsize; 143 void (*f)(void *ctx, int *nq, 144 const CeedScalar *u,const CeedScalar *u1, 145 const CeedScalar *u2,const CeedScalar *u3, 146 const CeedScalar *u4,const CeedScalar *u5, 147 const CeedScalar *u6,const CeedScalar *u7, 148 const CeedScalar *u8,const CeedScalar *u9, 149 const CeedScalar *u10,const CeedScalar *u11, 150 const CeedScalar *u12,const CeedScalar *u13, 151 const CeedScalar *u14,const CeedScalar *u15, 152 CeedScalar *v,CeedScalar *v1,CeedScalar *v2, 153 CeedScalar *v3,CeedScalar *v4,CeedScalar *v5, 154 CeedScalar *v6,CeedScalar *v7,CeedScalar *v8, 155 CeedScalar *v9, CeedScalar *v10,CeedScalar *v11, 156 CeedScalar *v12,CeedScalar *v13,CeedScalar *v14, 157 CeedScalar *v15, int *err); 158 } fContext; 159 160 struct CeedOperatorField_private { 161 CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity 162 CeedTransposeMode lmode; /// Transpose mode for lvector ordering 163 CeedBasis basis; /// Basis or NULL for collocated fields 164 CeedVector 165 vec; /// State vector for passive fields, NULL for active fields 166 }; 167 168 struct CeedOperator_private { 169 Ceed ceed; 170 int refcount; 171 int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 172 int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector, 173 CeedVector, CeedRequest *); 174 int (*Destroy)(CeedOperator); 175 CeedOperatorField *inputfields; 176 CeedOperatorField *outputfields; 177 CeedInt numelements; /// Number of elements 178 CeedInt numqpoints; /// Number of quadrature points over all elements 179 CeedInt nfields; /// Number of fields that have been set 180 CeedQFunction qf; 181 CeedQFunction dqf; 182 CeedQFunction dqfT; 183 bool setupdone; 184 void *data; 185 }; 186 187 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int, 188 const char *, va_list); 189 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int, 190 const char *, va_list); 191 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int, 192 const char *, va_list); 193 CEED_INTERN int CeedSetErrorHandler(Ceed ceed, 194 int (eh)(Ceed, const char *, int, const char *, 195 int, const char *, va_list)); 196 197 #endif 198