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