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 30 struct Ceed_private { 31 Ceed delegate; 32 int (*Error)(Ceed, const char *, int, const char *, int, const char *, va_list); 33 int (*Destroy)(Ceed); 34 int (*VecCreate)(CeedInt, CeedVector); 35 int (*ElemRestrictionCreate)(CeedMemType, CeedCopyMode, 36 const CeedInt *, CeedElemRestriction); 37 int (*ElemRestrictionCreateBlocked)(CeedMemType, CeedCopyMode, 38 const CeedInt *, CeedElemRestriction); 39 int (*BasisCreateTensorH1)(CeedInt, CeedInt, CeedInt, const CeedScalar *, 40 const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis); 41 int (*BasisCreateH1)(CeedElemTopology, CeedInt, CeedInt, CeedInt, 42 const CeedScalar *, 43 const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis); 44 int (*QFunctionCreate)(CeedQFunction); 45 int (*OperatorCreate)(CeedOperator); 46 int refcount; 47 void *data; 48 }; 49 50 struct CeedVector_private { 51 Ceed ceed; 52 int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *); 53 int (*SetValue)(CeedVector, CeedScalar); 54 int (*GetArray)(CeedVector, CeedMemType, CeedScalar **); 55 int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **); 56 int (*RestoreArray)(CeedVector, CeedScalar **); 57 int (*RestoreArrayRead)(CeedVector, const CeedScalar **); 58 int (*Destroy)(CeedVector); 59 int refcount; 60 CeedInt length; 61 uint64_t state; 62 void *data; 63 }; 64 65 struct CeedElemRestriction_private { 66 Ceed ceed; 67 int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedTransposeMode, 68 CeedVector, CeedVector, CeedRequest *); 69 int (*Destroy)(CeedElemRestriction); 70 int refcount; 71 CeedInt nelem; /* number of elements */ 72 CeedInt elemsize; /* number of dofs per element */ 73 CeedInt ndof; /* size of the L-vector, can be used for checking for 74 correct vector sizes */ 75 CeedInt ncomp; /* number of components */ 76 CeedInt blksize; /* number of elements in a batch */ 77 CeedInt nblk; /* number of blocks of elements */ 78 void *data; /* place for the backend to store any data */ 79 }; 80 81 struct CeedBasis_private { 82 Ceed ceed; 83 int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, 84 const CeedScalar *, 85 CeedScalar *); 86 int (*Destroy)(CeedBasis); 87 int refcount; 88 bool tensorbasis; /* flag for tensor basis */ 89 CeedInt dim; /* topological dimension */ 90 CeedInt ncomp; /* number of field components (1 for scalar fields) */ 91 CeedInt P1d; /* number of nodes in one dimension */ 92 CeedInt Q1d; /* number of quadrature points in one dimension */ 93 CeedInt P; /* total number of nodes */ 94 CeedInt Q; /* total number of quadrature points */ 95 CeedScalar *qref1d; /* Array of length Q1d holding the locations of 96 quadrature points on the 1D reference element [-1, 1] */ 97 CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on 98 the reference element */ 99 CeedScalar 100 *interp1d; /* row-major matrix of shape [Q1d, P1d] expressing the values of 101 nodal basis functions at quadrature points */ 102 CeedScalar 103 *grad1d; /* row-major matrix of shape [Q1d, P1d] matrix expressing derivatives of 104 nodal basis functions at quadrature points */ 105 void *data; /* place for the backend to store any data */ 106 }; 107 108 struct CeedQFunctionField_private { 109 const char *fieldname; 110 CeedInt ncomp; 111 CeedEvalMode emode; 112 }; 113 114 struct CeedQFunction_private { 115 Ceed ceed; 116 int (*Apply)(CeedQFunction, CeedInt, const CeedScalar *const *, 117 CeedScalar *const *); 118 int (*Destroy)(CeedQFunction); 119 int refcount; 120 CeedInt vlength; // Number of quadrature points must be padded to a multiple of vlength 121 CeedQFunctionField inputfields[16]; 122 CeedQFunctionField outputfields[16]; 123 CeedInt numinputfields, numoutputfields; 124 int (*function)(void*, CeedInt, const CeedScalar *const*, CeedScalar *const*); 125 const char *focca; 126 void *ctx; /* user context for function */ 127 size_t ctxsize; /* size of user context; may be used to copy to a device */ 128 void *data; /* backend data */ 129 }; 130 131 struct CeedOperatorField_private { 132 CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity 133 CeedBasis basis; /// Basis or NULL for collocated fields 134 CeedVector 135 vec; /// State vector for passive fields, NULL for active fields 136 }; 137 138 struct CeedOperator_private { 139 Ceed ceed; 140 int refcount; 141 int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 142 int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector, 143 CeedVector, CeedRequest *); 144 int (*GetQData)(CeedOperator, CeedVector *); 145 int (*Destroy)(CeedOperator); 146 CeedOperatorField inputfields[16]; 147 CeedOperatorField outputfields[16]; 148 CeedInt numelements; /// Number of elements 149 CeedInt numqpoints; /// Number of quadrature points over all elements 150 CeedInt nfields; /// Number of fields that have been set 151 CeedQFunction qf; 152 CeedQFunction dqf; 153 CeedQFunction dqfT; 154 bool setupdone; 155 void *data; 156 }; 157 158 CEED_INTERN int CeedErrorReturn(Ceed, const char *, int, const char *, int, 159 const char *, va_list); 160 CEED_INTERN int CeedErrorAbort(Ceed, const char *, int, const char *, int, 161 const char *, va_list); 162 CEED_INTERN int CeedErrorExit(Ceed, const char *, int, const char *, int, 163 const char *, va_list); 164 CEED_INTERN int CeedSetErrorHandler(Ceed ceed, 165 int (eh)(Ceed, const char *, int, const char *, 166 int, const char *, va_list)); 167 168 #endif 169