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