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/ceed.h> 23 #include <ceed/backend.h> 24 #include <stdbool.h> 25 26 /** @defgroup CeedUser Public API for Ceed 27 @ingroup Ceed 28 */ 29 /** @defgroup CeedBackend Backend API for Ceed 30 @ingroup Ceed 31 */ 32 /** @defgroup CeedDeveloper Internal library functions for Ceed 33 @ingroup Ceed 34 */ 35 /** @defgroup CeedVectorUser Public API for CeedVector 36 @ingroup CeedVector 37 */ 38 /** @defgroup CeedVectorBackend Backend API for CeedVector 39 @ingroup CeedVector 40 */ 41 /** @defgroup CeedVectorDeveloper Internal library functions for CeedVector 42 @ingroup CeedVector 43 */ 44 /** @defgroup CeedElemRestrictionUser Public API for CeedElemRestriction 45 @ingroup CeedElemRestriction 46 */ 47 /** @defgroup CeedElemRestrictionBackend Backend API for CeedElemRestriction 48 @ingroup CeedElemRestriction 49 */ 50 /** @defgroup CeedElemRestrictionDeveloper Internal library functions for CeedElemRestriction 51 @ingroup CeedElemRestriction 52 */ 53 /** @defgroup CeedBasisUser Public API for CeedBasis 54 @ingroup CeedBasis 55 */ 56 /** @defgroup CeedBasisBackend Backend API for CeedBasis 57 @ingroup CeedBasis 58 */ 59 /** @defgroup CeedBasisDeveloper Internal library functions for CeedBasis 60 @ingroup CeedBasis 61 */ 62 /** @defgroup CeedQFunctionUser Public API for CeedQFunction 63 @ingroup CeedQFunction 64 */ 65 /** @defgroup CeedQFunctionBackend Backend API for CeedQFunction 66 @ingroup CeedQFunction 67 */ 68 /** @defgroup CeedQFunctionDeveloper Internal library functions for CeedQFunction 69 @ingroup CeedQFunction 70 */ 71 /** @defgroup CeedOperatorUser Public API for CeedOperator 72 @ingroup CeedOperator 73 */ 74 /** @defgroup CeedOperatorBackend Backend API for CeedOperator 75 @ingroup CeedOperator 76 */ 77 /** @defgroup CeedOperatorDeveloper Internal library functions for CeedOperator 78 @ingroup CeedOperator 79 */ 80 81 // Lookup table field for backend functions 82 typedef struct { 83 const char *func_name; 84 size_t offset; 85 } FOffset; 86 87 // Lookup table field for object delegates 88 typedef struct { 89 char *obj_name; 90 Ceed delegate; 91 } ObjDelegate; 92 93 struct Ceed_private { 94 const char *resource; 95 Ceed delegate; 96 Ceed parent; 97 ObjDelegate *obj_delegates; 98 int obj_delegate_count; 99 Ceed op_fallback_ceed, op_fallback_parent; 100 const char *op_fallback_resource; 101 int (*Error)(Ceed, const char *, int, const char *, int, const char *, 102 va_list *); 103 int (*GetPreferredMemType)(CeedMemType *); 104 int (*Destroy)(Ceed); 105 int (*VectorCreate)(CeedInt, CeedVector); 106 int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode, 107 const CeedInt *, CeedElemRestriction); 108 int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode, 109 const CeedInt *, CeedElemRestriction); 110 int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *, 111 const CeedScalar *, const CeedScalar *, 112 const CeedScalar *, CeedBasis); 113 int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt, 114 const CeedScalar *, 115 const CeedScalar *, const CeedScalar *, 116 const CeedScalar *, CeedBasis); 117 int (*TensorContractCreate)(CeedBasis, CeedTensorContract); 118 int (*QFunctionCreate)(CeedQFunction); 119 int (*QFunctionContextCreate)(CeedQFunctionContext); 120 int (*OperatorCreate)(CeedOperator); 121 int (*CompositeOperatorCreate)(CeedOperator); 122 int ref_count; 123 bool is_deterministic; 124 void *data; 125 bool debug; 126 char err_msg[CEED_MAX_RESOURCE_LEN]; 127 FOffset *f_offsets; 128 }; 129 130 struct CeedVector_private { 131 Ceed ceed; 132 int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *); 133 int (*SetValue)(CeedVector, CeedScalar); 134 int (*SyncArray)(CeedVector, CeedMemType); 135 int (*TakeArray)(CeedVector, CeedMemType, CeedScalar **); 136 int (*GetArray)(CeedVector, CeedMemType, CeedScalar **); 137 int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **); 138 int (*RestoreArray)(CeedVector); 139 int (*RestoreArrayRead)(CeedVector); 140 int (*Norm)(CeedVector, CeedNormType, CeedScalar *); 141 int (*Reciprocal)(CeedVector); 142 int (*Destroy)(CeedVector); 143 int ref_count; 144 CeedInt length; 145 uint64_t state; 146 uint64_t num_readers; 147 void *data; 148 }; 149 150 struct CeedElemRestriction_private { 151 Ceed ceed; 152 int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedVector, CeedVector, 153 CeedRequest *); 154 int (*ApplyBlock)(CeedElemRestriction, CeedInt, CeedTransposeMode, CeedVector, 155 CeedVector, CeedRequest *); 156 int (*GetOffsets)(CeedElemRestriction, CeedMemType, const CeedInt **); 157 int (*Destroy)(CeedElemRestriction); 158 int ref_count; 159 CeedInt num_elem; /* number of elements */ 160 CeedInt elem_size; /* number of nodes per element */ 161 CeedInt num_comp; /* number of components */ 162 CeedInt comp_stride; /* Component stride for L-vector ordering */ 163 CeedInt l_size; /* size of the L-vector, can be used for checking 164 for correct vector sizes */ 165 CeedInt blk_size; /* number of elements in a batch */ 166 CeedInt num_blk; /* number of blocks of elements */ 167 CeedInt *strides; /* strides between [nodes, components, elements] */ 168 CeedInt layout[3]; /* E-vector layout [nodes, components, elements] */ 169 uint64_t num_readers; /* number of instances of offset read only access */ 170 void *data; /* place for the backend to store any data */ 171 }; 172 173 struct CeedBasis_private { 174 Ceed ceed; 175 int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, 176 CeedVector, CeedVector); 177 int (*Destroy)(CeedBasis); 178 int ref_count; 179 bool tensor_basis; /* flag for tensor basis */ 180 CeedInt dim; /* topological dimension */ 181 CeedElemTopology topo; /* element topology */ 182 CeedInt num_comp; /* number of field components (1 for scalar fields) */ 183 CeedInt P_1d; /* number of nodes in one dimension */ 184 CeedInt Q_1d; /* number of quadrature points in one dimension */ 185 CeedInt P; /* total number of nodes */ 186 CeedInt Q; /* total number of quadrature points */ 187 CeedScalar *q_ref_1d; /* Array of length Q1d holding the locations of 188 quadrature points on the 1D reference 189 element [-1, 1] */ 190 CeedScalar 191 *q_weight_1d; /* array of length Q1d holding the quadrature weights on 192 the reference element */ 193 CeedScalar 194 *interp; /* row-major matrix of shape [Q, P] expressing the values of 195 nodal basis functions at quadrature points */ 196 CeedScalar 197 *interp_1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of 198 nodal basis functions at quadrature points */ 199 CeedScalar 200 *grad; /* row-major matrix of shape [dim*Q, P] matrix expressing 201 derivatives of nodal basis functions at quadrature points */ 202 CeedScalar 203 *grad_1d; /* row-major matrix of shape [Q1d, P1d] matrix expressing 204 derivatives of nodal basis functions at quadrature points */ 205 CeedTensorContract contract; /* tensor contraction object */ 206 void *data; /* place for the backend to store any data */ 207 }; 208 209 struct CeedTensorContract_private { 210 Ceed ceed; 211 int (*Apply)(CeedTensorContract, CeedInt, CeedInt, CeedInt, CeedInt, 212 const CeedScalar *restrict, CeedTransposeMode, const CeedInt, 213 const CeedScalar *restrict, CeedScalar *restrict); 214 int (*Destroy)(CeedTensorContract); 215 int ref_count; 216 void *data; 217 }; 218 219 struct CeedQFunctionField_private { 220 const char *field_name; 221 CeedInt size; 222 CeedEvalMode eval_mode; 223 }; 224 225 struct CeedQFunction_private { 226 Ceed ceed; 227 int (*Apply)(CeedQFunction, CeedInt, CeedVector *, CeedVector *); 228 int (*SetCUDAUserFunction)(CeedQFunction, void *); 229 int (*SetHIPUserFunction)(CeedQFunction, void *); 230 int (*Destroy)(CeedQFunction); 231 int ref_count; 232 CeedInt vec_length; /* Number of quadrature points must be padded to a 233 multiple of vec_length */ 234 CeedQFunctionField *input_fields; 235 CeedQFunctionField *output_fields; 236 CeedInt num_input_fields, num_output_fields; 237 CeedQFunctionUser function; 238 const char *source_path; 239 const char *qf_name; 240 bool identity; 241 bool fortran_status; 242 CeedInt operators_set; 243 CeedQFunctionContext ctx; /* user context for function */ 244 void *data; /* place for the backend to store any data */ 245 }; 246 247 struct CeedQFunctionContext_private { 248 Ceed ceed; 249 int ref_count; 250 int (*SetData)(CeedQFunctionContext, CeedMemType, CeedCopyMode, void *); 251 int (*GetData)(CeedQFunctionContext, CeedMemType, void *); 252 int (*RestoreData)(CeedQFunctionContext); 253 int (*Destroy)(CeedQFunctionContext); 254 uint64_t state; 255 size_t ctx_size; 256 void *data; 257 }; 258 259 /// Struct to handle the context data to use the Fortran QFunction stub 260 /// @ingroup CeedQFunction 261 struct CeedFortranContext_private { 262 CeedQFunctionContext innerctx; 263 void (*f)(void *ctx, int *nq, 264 const CeedScalar *u,const CeedScalar *u1, 265 const CeedScalar *u2,const CeedScalar *u3, 266 const CeedScalar *u4,const CeedScalar *u5, 267 const CeedScalar *u6,const CeedScalar *u7, 268 const CeedScalar *u8,const CeedScalar *u9, 269 const CeedScalar *u10,const CeedScalar *u11, 270 const CeedScalar *u12,const CeedScalar *u13, 271 const CeedScalar *u14,const CeedScalar *u15, 272 CeedScalar *v,CeedScalar *v1,CeedScalar *v2, 273 CeedScalar *v3,CeedScalar *v4,CeedScalar *v5, 274 CeedScalar *v6,CeedScalar *v7,CeedScalar *v8, 275 CeedScalar *v9, CeedScalar *v10,CeedScalar *v11, 276 CeedScalar *v12,CeedScalar *v13,CeedScalar *v14, 277 CeedScalar *v15, int *err); 278 }; 279 typedef struct CeedFortranContext_private *CeedFortranContext; 280 281 struct CeedOperatorField_private { 282 CeedElemRestriction elem_restr; /* Restriction from L-vector */ 283 CeedBasis basis; /* Basis or CEED_BASIS_COLLOCATED for 284 collocated fields */ 285 CeedVector vec; /* State vector for passive fields or 286 CEED_VECTOR_NONE for no vector */ 287 const char *field_name; /* matching QFunction field name */ 288 }; 289 290 struct CeedOperator_private { 291 Ceed ceed; 292 CeedOperator op_fallback; 293 CeedQFunction qf_fallback; 294 int ref_count; 295 int (*LinearAssembleQFunction)(CeedOperator, CeedVector *, 296 CeedElemRestriction *, CeedRequest *); 297 int (*LinearAssembleDiagonal)(CeedOperator, CeedVector, CeedRequest *); 298 int (*LinearAssembleAddDiagonal)(CeedOperator, CeedVector, CeedRequest *); 299 int (*LinearAssemblePointBlockDiagonal)(CeedOperator, CeedVector, 300 CeedRequest *); 301 int (*LinearAssembleAddPointBlockDiagonal)(CeedOperator, CeedVector, 302 CeedRequest *); 303 int (*LinearAssembleSymbolic)(CeedOperator, CeedInt *, CeedInt **, CeedInt **); 304 int (*LinearAssemble)(CeedOperator, CeedVector); 305 int (*CreateFDMElementInverse)(CeedOperator, CeedOperator *, CeedRequest *); 306 int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 307 int (*ApplyComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 308 int (*ApplyAdd)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 309 int (*ApplyAddComposite)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 310 int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector, 311 CeedVector, CeedRequest *); 312 int (*Destroy)(CeedOperator); 313 CeedOperatorField *input_fields; 314 CeedOperatorField *output_fields; 315 CeedInt num_elem; /* Number of elements */ 316 CeedInt num_qpts; /* Number of quadrature points over all elements */ 317 CeedInt num_fields; /* Number of fields that have been set */ 318 CeedQFunction qf; 319 CeedQFunction dqf; 320 CeedQFunction dqfT; 321 bool interface_setup; 322 bool backend_setup; 323 bool composite; 324 bool has_restriction; 325 CeedOperator *sub_operators; 326 CeedInt num_suboperators; 327 void *data; 328 }; 329 330 #endif 331