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 #ifndef _ceed_impl_h 18 #define _ceed_impl_h 19 20 #include <ceed.h> 21 #include <stdbool.h> 22 23 #define CEED_INTERN CEED_EXTERN __attribute__((visibility ("hidden"))) 24 25 #define CEED_MAX_RESOURCE_LEN 1024 26 #define CEED_ALIGN 64 27 28 struct Ceed_private { 29 int (*Error)(Ceed, const char *, int, const char *, int, const char *, va_list); 30 int (*Destroy)(Ceed); 31 int (*VecCreate)(Ceed, CeedInt, CeedVector); 32 int (*ElemRestrictionCreate)(CeedElemRestriction, CeedMemType, CeedCopyMode, 33 const CeedInt *); 34 int (*ElemRestrictionCreateBlocked)(CeedElemRestriction, CeedMemType, CeedCopyMode, 35 const CeedInt *); 36 int (*BasisCreateTensorH1)(Ceed, CeedInt, CeedInt, CeedInt, const CeedScalar *, 37 const CeedScalar *, const CeedScalar *, const CeedScalar *, CeedBasis); 38 int (*QFunctionCreate)(CeedQFunction); 39 int (*OperatorCreate)(CeedOperator); 40 int refcount; 41 void *data; 42 }; 43 44 /* In the next 3 functions, p has to be the address of a pointer type, i.e. p 45 has to be a pointer to a pointer. */ 46 CEED_INTERN int CeedMallocArray(size_t n, size_t unit, void *p); 47 CEED_INTERN int CeedCallocArray(size_t n, size_t unit, void *p); 48 CEED_INTERN int CeedReallocArray(size_t n, size_t unit, void *p); 49 CEED_INTERN int CeedFree(void *p); 50 51 #define CeedChk(ierr) do { if (ierr) return ierr; } while (0) 52 /* Note that CeedMalloc and CeedCalloc will, generally, return pointers with 53 different memory alignments: CeedMalloc returns pointers aligned at 54 CEED_ALIGN bytes, while CeedCalloc uses the alignment of calloc. */ 55 #define CeedMalloc(n, p) CeedMallocArray((n), sizeof(**(p)), p) 56 #define CeedCalloc(n, p) CeedCallocArray((n), sizeof(**(p)), p) 57 #define CeedRealloc(n, p) CeedReallocArray((n), sizeof(**(p)), p) 58 59 struct CeedVector_private { 60 Ceed ceed; 61 int (*SetArray)(CeedVector, CeedMemType, CeedCopyMode, CeedScalar *); 62 int (*SetValue)(CeedVector, CeedScalar); 63 int (*GetArray)(CeedVector, CeedMemType, CeedScalar **); 64 int (*GetArrayRead)(CeedVector, CeedMemType, const CeedScalar **); 65 int (*RestoreArray)(CeedVector, CeedScalar **); 66 int (*RestoreArrayRead)(CeedVector, const CeedScalar **); 67 int (*Destroy)(CeedVector); 68 int refcount; 69 CeedInt length; 70 void *data; 71 }; 72 73 struct CeedElemRestriction_private { 74 Ceed ceed; 75 int (*Apply)(CeedElemRestriction, CeedTransposeMode, CeedTransposeMode, 76 CeedVector, CeedVector, CeedRequest *); 77 int (*Destroy)(CeedElemRestriction); 78 int refcount; 79 CeedInt nelem; /* number of elements */ 80 CeedInt elemsize; /* number of dofs per element */ 81 CeedInt ndof; /* size of the L-vector, can be used for checking for 82 correct vector sizes */ 83 CeedInt ncomp; /* number of components */ 84 CeedInt blksize; /* number of elements in a batch */ 85 CeedInt nblk; /* number of blocks of elements */ 86 void *data; /* place for the backend to store any data */ 87 }; 88 89 struct CeedBasis_private { 90 Ceed ceed; 91 int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, const CeedScalar *, 92 CeedScalar *); 93 int (*Destroy)(CeedBasis); 94 int refcount; 95 CeedInt dim; /* topological dimension */ 96 CeedInt ncomp; /* number of field components (1 for scalar fields) */ 97 CeedInt P1d; /* number of nodes in one dimension */ 98 CeedInt Q1d; /* number of quadrature points in one dimension */ 99 CeedScalar *qref1d; /* Array of length Q1d holding the locations of 100 quadrature points on the 1D reference element [-1, 1] */ 101 CeedScalar *qweight1d; /* array of length Q1d holding the quadrature weights on 102 the reference element */ 103 CeedScalar *interp1d; /* row-major Q1d × P1d matrix expressing the values of 104 nodal basis functions at quadrature points */ 105 CeedScalar *grad1d; /* row-major Q1d × P1d matrix expressing derivatives of 106 nodal basis functions at quadrature points */ 107 void *data; /* place for the backend to store any data */ 108 }; 109 110 struct CeedQFunctionField { 111 const char *fieldname; 112 CeedInt ncomp; 113 CeedEvalMode emode; 114 }; 115 116 struct CeedQFunction_private { 117 Ceed ceed; 118 int (*Apply)(CeedQFunction, CeedInt, const CeedScalar *const *, 119 CeedScalar *const *); 120 int (*Destroy)(CeedQFunction); 121 int refcount; 122 CeedInt vlength; // Number of quadrature points must be padded to a multiple of vlength 123 struct CeedQFunctionField inputfields[16]; 124 struct CeedQFunctionField outputfields[16]; 125 CeedInt numinputfields, numoutputfields; 126 int (*function)(void*, CeedInt, const CeedScalar *const*, CeedScalar *const*); 127 const char *focca; 128 void *ctx; /* user context for function */ 129 size_t ctxsize; /* size of user context; may be used to copy to a device */ 130 void *data; /* backend data */ 131 }; 132 133 struct CeedOperatorField { 134 CeedElemRestriction Erestrict; /// Restriction from L-vector or NULL if identity 135 CeedBasis basis; /// Basis or NULL for collocated fields 136 CeedVector vec; /// State vector for passive fields, NULL for active fields 137 }; 138 139 struct CeedOperator_private { 140 Ceed ceed; 141 int refcount; 142 int (*Apply)(CeedOperator, CeedVector, CeedVector, CeedRequest *); 143 int (*ApplyJacobian)(CeedOperator, CeedVector, CeedVector, CeedVector, 144 CeedVector, CeedRequest *); 145 int (*GetQData)(CeedOperator, CeedVector *); 146 int (*Destroy)(CeedOperator); 147 struct CeedOperatorField inputfields[16]; 148 struct CeedOperatorField outputfields[16]; 149 CeedInt numelements; /// Number of elements 150 CeedInt numqpoints; /// Number of quadrature points over all elements 151 CeedInt nfields; /// Number of fields that have been set 152 CeedQFunction qf; 153 CeedQFunction dqf; 154 CeedQFunction dqfT; 155 bool setupdone; 156 void *data; 157 }; 158 159 static inline CeedInt CeedIntMin(CeedInt a, CeedInt b) { return a < b ? a : b; } 160 161 #endif 162