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, CeedTransposeMode, 92 CeedVector, CeedVector, CeedRequest *); 93 int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, 94 CeedTransposeMode, CeedVector, CeedVector, CeedRequest *); 95 int (*Destroy)(CeedElemRestriction); 96 int refcount; 97 CeedInt nelem; /* number of elements */ 98 CeedInt elemsize; /* number of nodes per element */ 99 CeedInt nnodes; /* size of the L-vector, can be used for checking for 100 correct vector sizes */ 101 CeedInt ncomp; /* number of components */ 102 CeedInt blksize; /* number of elements in a batch */ 103 CeedInt nblk; /* number of blocks of elements */ 104 void *data; /* place for the backend to store any data */ 105 }; 106 107 struct CeedBasis_private { 108 Ceed ceed; 109 int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, 110 CeedVector, CeedVector); 111 int (*Destroy)(CeedBasis); 112 int refcount; 113 bool tensorbasis; /* flag for tensor basis */ 114 CeedInt dim; /* topological dimension */ 115 CeedInt ncomp; /* number of field components (1 for scalar fields) */ 116 CeedInt P1d; /* number of nodes in one dimension */ 117 CeedInt Q1d; /* number of quadrature points in one dimension */ 118 CeedInt P; /* total number of nodes */ 119 CeedInt Q; /* total number of quadrature points */ 120 CeedScalar *qref1d; /* Array of length Q1d holding the locations of 121 quadrature points on the 1D reference element [-1, 1] */ 122 CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on 123 the reference element */ 124 CeedScalar 125 *interp; /* row-major matrix of shape [Q, P] expressing the values of 126 nodal basis functions at quadrature points */ 127 CeedScalar 128 *interp1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of 129 nodal basis functions at quadrature points */ 130 CeedScalar 131 *grad; /* row-major matrix of shape [dim*Q, P] matrix expressing derivatives of 132 nodal basis functions at quadrature points */ 133 CeedScalar 134 *grad1d; /* row-major matrix of shape [Q1d, P1d] matrix expressing derivatives of 135 nodal basis functions at quadrature points */ 136 CeedTensorContract contract; /* tensor contraction object */ 137 void *data; /* place for the backend to store any data */ 138 }; 139 140 struct CeedTensorContract_private { 141 Ceed ceed; 142 int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt, 143 const CeedScalar *restrict, CeedTransposeMode, const CeedInt, 144 const CeedScalar *restrict, CeedScalar *restrict); 145 int (*Destroy)(CeedTensorContract); 146 int refcount; 147 void *data; 148 }; 149 150 struct CeedQFunctionField_private { 151 const char *fieldname; 152 CeedInt size; 153 CeedEvalMode emode; 154 }; 155 156 struct CeedQFunction_private { 157 Ceed ceed; 158 int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *); 159 int (*Destroy)(CeedQFunction); 160 int refcount; 161 CeedInt vlength; // Number of quadrature points must be padded to a multiple of vlength 162 CeedQFunctionField *inputfields; 163 CeedQFunctionField *outputfields; 164 CeedInt numinputfields, numoutputfields; 165 CeedQFunctionUser function; 166 const char *sourcepath; 167 const char *qfname; 168 bool fortranstatus; 169 bool identity; 170 void *ctx; /* user context for function */ 171 size_t ctxsize; /* size of user context; may be used to copy to a device */ 172 void *data; /* backend data */ 173 }; 174 175 /// Struct to handle the context data to use the Fortran QFunction stub 176 /// @ingroup CeedQFunction 177 typedef struct { 178 CeedScalar *innerctx; 179 size_t innerctxsize; 180 void (*f)(void *ctx, int *nq, 181 const CeedScalar *u,const CeedScalar *u1, 182 const CeedScalar *u2,const CeedScalar *u3, 183 const CeedScalar *u4,const CeedScalar *u5, 184 const CeedScalar *u6,const CeedScalar *u7, 185 const CeedScalar *u8,const CeedScalar *u9, 186 const CeedScalar *u10,const CeedScalar *u11, 187 const CeedScalar *u12,const CeedScalar *u13, 188 const CeedScalar *u14,const CeedScalar *u15, 189 CeedScalar *v,CeedScalar *v1,CeedScalar *v2, 190 CeedScalar *v3,CeedScalar *v4,CeedScalar *v5, 191 CeedScalar *v6,CeedScalar *v7,CeedScalar *v8, 192 CeedScalar *v9, CeedScalar *v10,CeedScalar *v11, 193 CeedScalar *v12,CeedScalar *v13,CeedScalar *v14, 194 CeedScalar *v15, int *err); 195 } fContext; 196 197 struct CeedOperatorField_private { 198 CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity 199 CeedTransposeMode lmode; /// Transpose mode for lvector ordering 200 CeedBasis basis; /// Basis or NULL for collocated fields 201 CeedVector 202 vec; /// State vector for passive fields, NULL for active fields 203 }; 204 205 struct CeedOperator_private { 206 Ceed ceed; 207 CeedOperator opfallback; 208 CeedQFunction qffallback; 209 int refcount; 210 int (*AssembleLinearQFunction)(CeedOperator, CeedVector *, 211 CeedElemRestriction *, CeedRequest *); 212 int (*AssembleLinearDiagonal)(CeedOperator, CeedVector *, CeedRequest *); 213 int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *); 214 int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 215 int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 216 int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 217 int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 218 int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector, 219 CeedVector, CeedRequest *); 220 int (*Destroy)(CeedOperator); 221 CeedOperatorField *inputfields; 222 CeedOperatorField *outputfields; 223 CeedInt numelements; /// Number of elements 224 CeedInt numqpoints; /// Number of quadrature points over all elements 225 CeedInt nfields; /// Number of fields that have been set 226 CeedQFunction qf; 227 CeedQFunction dqf; 228 CeedQFunction dqfT; 229 bool setupdone; 230 bool composite; 231 bool hasrestriction; 232 CeedOperator *suboperators; 233 CeedInt numsub; 234 void *data; 235 }; 236 237 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int, 238 const char *, va_list); 239 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int, 240 const char *, va_list); 241 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int, 242 const char *, va_list); 243 CEED_INTERN int CeedSetErrorHandler(Ceed ceed, 244 int (eh)(Ceed, const char *, int, 245 const char *, int, const char *, 246 va_list)); 247 248 #endif 249