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 <ceed-backend.h> 24 #include <stdbool.h> 25 26 // Lookup table field for backend functions 27 typedef struct { 28 const char *fname; 29 size_t offset; 30 } foffset; 31 32 // Lookup table field for object delegates 33 typedef struct { 34 char *objname; 35 Ceed delegate; 36 } objdelegate; 37 38 struct Ceed_private { 39 const char *resource; 40 Ceed delegate; 41 Ceed parent; 42 objdelegate *objdelegates; 43 int objdelegatecount; 44 Ceed opfallbackceed, opfallbackparent; 45 const char *opfallbackresource; 46 int (*Error)(Ceed, const char *, int, const char *, int, const char *, 47 va_list); 48 int (*GetPreferredMemType)(CeedMemType *); 49 int (*Destroy)(Ceed); 50 int (*VectorCreate)(CeedInt, CeedVector); 51 int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode, 52 const CeedInt *, CeedElemRestriction); 53 int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode, 54 const CeedInt *, CeedElemRestriction); 55 int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *, 56 const CeedScalar *, const CeedScalar *, 57 const CeedScalar *, CeedBasis); 58 int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt, 59 const CeedScalar *, 60 const CeedScalar *, const CeedScalar *, 61 const CeedScalar *, CeedBasis); 62 int (*TensorContractCreate)(CeedBasis, CeedTensorContract); 63 int (*QFunctionCreate)(CeedQFunction); 64 int (*OperatorCreate)(CeedOperator); 65 int (*CompositeOperatorCreate)(CeedOperator); 66 int refcount; 67 void *data; 68 foffset *foffsets; 69 }; 70 71 struct CeedVector_private { 72 Ceed ceed; 73 int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *); 74 int (*SetValue)(CeedVector, CeedScalar); 75 int (*SyncArray)(CeedVector, CeedMemType); 76 int (*GetArray)(CeedVector, CeedMemType, CeedScalar **); 77 int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **); 78 int (*RestoreArray)(CeedVector); 79 int (*RestoreArrayRead)(CeedVector); 80 int (*Norm)(CeedVector, CeedNormType, CeedScalar *); 81 int (*Destroy)(CeedVector); 82 int refcount; 83 CeedInt length; 84 uint64_t state; 85 uint64_t numreaders; 86 void *data; 87 }; 88 89 struct CeedElemRestriction_private { 90 Ceed ceed; 91 int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector, 92 CeedRequest *); 93 int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector, 94 CeedVector, CeedRequest *); 95 int (*Destroy)(CeedElemRestriction); 96 int refcount; 97 CeedTransposeMode lmode; /* Transpose mode for L-vector ordering */ 98 CeedInt nelem; /* number of elements */ 99 CeedInt elemsize; /* number of nodes per element */ 100 CeedInt nnodes; /* size of the L-vector, can be used for checking 101 for correct vector sizes */ 102 CeedInt ncomp; /* number of components */ 103 CeedInt blksize; /* number of elements in a batch */ 104 CeedInt nblk; /* number of blocks of elements */ 105 void *data; /* place for the backend to store any data */ 106 }; 107 108 struct CeedBasis_private { 109 Ceed ceed; 110 int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, 111 CeedVector, CeedVector); 112 int (*Destroy)(CeedBasis); 113 int refcount; 114 bool tensorbasis; /* flag for tensor basis */ 115 CeedInt dim; /* topological dimension */ 116 CeedInt ncomp; /* number of field components (1 for scalar fields) */ 117 CeedInt P1d; /* number of nodes in one dimension */ 118 CeedInt Q1d; /* number of quadrature points in one dimension */ 119 CeedInt P; /* total number of nodes */ 120 CeedInt Q; /* total number of quadrature points */ 121 CeedScalar *qref1d; /* Array of length Q1d holding the locations of 122 quadrature points on the 1D reference 123 element [-1, 1] */ 124 CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on 125 the reference element */ 126 CeedScalar 127 *interp; /* row-major matrix of shape [Q, P] expressing the values of 128 nodal basis functions at quadrature points */ 129 CeedScalar 130 *interp1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of 131 nodal basis functions at quadrature points */ 132 CeedScalar 133 *grad; /* row-major matrix of shape [dim*Q, P] matrix expressing 134 derivatives of nodal basis functions at quadrature points */ 135 CeedScalar 136 *grad1d; /* row-major matrix of shape [Q1d, P1d] matrix expressing 137 derivatives of nodal basis functions at quadrature points */ 138 CeedTensorContract contract; /* tensor contraction object */ 139 void *data; /* place for the backend to store any data */ 140 }; 141 142 struct CeedTensorContract_private { 143 Ceed ceed; 144 int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt, 145 const CeedScalar *restrict, CeedTransposeMode, const CeedInt, 146 const CeedScalar *restrict, CeedScalar *restrict); 147 int (*Destroy)(CeedTensorContract); 148 int refcount; 149 void *data; 150 }; 151 152 struct CeedQFunctionField_private { 153 const char *fieldname; 154 CeedInt size; 155 CeedEvalMode emode; 156 }; 157 158 struct CeedQFunction_private { 159 Ceed ceed; 160 int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *); 161 int (*Destroy)(CeedQFunction); 162 int refcount; 163 CeedInt vlength; /* Number of quadrature points must be padded to a 164 multiple of vlength */ 165 CeedQFunctionField *inputfields; 166 CeedQFunctionField *outputfields; 167 CeedInt numinputfields, numoutputfields; 168 CeedQFunctionUser function; 169 const char *sourcepath; 170 const char *qfname; 171 bool fortranstatus; 172 bool identity; 173 void *ctx; /* user context for function */ 174 size_t ctxsize; /* size of user context; may be used to copy to a device */ 175 void *data; /* place for the backend to store any data */ 176 }; 177 178 /// Struct to handle the context data to use the Fortran QFunction stub 179 /// @ingroup CeedQFunction 180 typedef struct { 181 CeedScalar *innerctx; 182 size_t innerctxsize; 183 void (*f)(void *ctx, int *nq, 184 const CeedScalar *u,const CeedScalar *u1, 185 const CeedScalar *u2,const CeedScalar *u3, 186 const CeedScalar *u4,const CeedScalar *u5, 187 const CeedScalar *u6,const CeedScalar *u7, 188 const CeedScalar *u8,const CeedScalar *u9, 189 const CeedScalar *u10,const CeedScalar *u11, 190 const CeedScalar *u12,const CeedScalar *u13, 191 const CeedScalar *u14,const CeedScalar *u15, 192 CeedScalar *v,CeedScalar *v1,CeedScalar *v2, 193 CeedScalar *v3,CeedScalar *v4,CeedScalar *v5, 194 CeedScalar *v6,CeedScalar *v7,CeedScalar *v8, 195 CeedScalar *v9, CeedScalar *v10,CeedScalar *v11, 196 CeedScalar *v12,CeedScalar *v13,CeedScalar *v14, 197 CeedScalar *v15, int *err); 198 } fContext; 199 200 struct CeedOperatorField_private { 201 CeedElemRestriction Erestrict; /* Restriction from L-vector */ 202 CeedBasis basis; /* Basis or CEED_BASIS_COLLOCATED for 203 collocated fields */ 204 CeedVector vec; /* State vector for passive fields or 205 CEED_VECTOR_NONE for no vector */ 206 }; 207 208 struct CeedOperator_private { 209 Ceed ceed; 210 CeedOperator opfallback; 211 CeedQFunction qffallback; 212 int refcount; 213 int (*AssembleLinearQFunction)(CeedOperator, CeedVector *, 214 CeedElemRestriction *, CeedRequest *); 215 int (*AssembleLinearDiagonal)(CeedOperator, CeedVector *, CeedRequest *); 216 int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *); 217 int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 218 int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 219 int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 220 int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 221 int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector, 222 CeedVector, CeedRequest *); 223 int (*Destroy)(CeedOperator); 224 CeedOperatorField *inputfields; 225 CeedOperatorField *outputfields; 226 CeedInt numelements; /// Number of elements 227 CeedInt numqpoints; /// Number of quadrature points over all elements 228 CeedInt nfields; /// Number of fields that have been set 229 CeedQFunction qf; 230 CeedQFunction dqf; 231 CeedQFunction dqfT; 232 bool setupdone; 233 bool composite; 234 bool hasrestriction; 235 CeedOperator *suboperators; 236 CeedInt numsub; 237 void *data; 238 }; 239 240 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int, 241 const char *, va_list); 242 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int, 243 const char *, va_list); 244 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int, 245 const char *, va_list); 246 CEED_INTERN int CeedSetErrorHandler(Ceed ceed, 247 int (eh)(Ceed, const char *, int, 248 const char *, int, const char *, 249 va_list)); 250 251 #endif 252