1 // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC. 2 // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707. 3 // All Rights 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 #include <ceed-impl.h> 18 #include <string.h> 19 #include "ceed-ref.h" 20 21 // Contracts on the middle index 22 // NOTRANSPOSE: V_ajc = T_jb U_abc 23 // TRANSPOSE: V_ajc = T_bj U_abc 24 // If Add != 0, "=" is replaced by "+=" 25 static int CeedTensorContract_Ref(Ceed ceed, 26 CeedInt A, CeedInt B, CeedInt C, CeedInt J, 27 const CeedScalar *restrict t, CeedTransposeMode tmode, 28 const CeedInt Add, 29 const CeedScalar *restrict u, CeedScalar *restrict v) { 30 CeedInt tstride0 = B, tstride1 = 1; 31 if (tmode == CEED_TRANSPOSE) { 32 tstride0 = 1; tstride1 = J; 33 } 34 35 if (!Add) 36 for (CeedInt q=0; q<A*J*C; q++) 37 v[q] = (CeedScalar) 0.0; 38 39 for (CeedInt a=0; a<A; a++) 40 for (CeedInt b=0; b<B; b++) 41 for (CeedInt j=0; j<J; j++) { 42 CeedScalar tq = t[j*tstride0 + b*tstride1]; 43 for (CeedInt c=0; c<C; c++) 44 v[(a*J+j)*C+c] += tq * u[(a*B+b)*C+c]; 45 } 46 return 0; 47 } 48 49 static int CeedBasisApply_Ref(CeedBasis basis, CeedInt nelem, 50 CeedTransposeMode tmode, CeedEvalMode emode, 51 const CeedScalar *u, CeedScalar *v) { 52 int ierr; 53 const CeedInt dim = basis->dim; 54 const CeedInt ncomp = basis->ncomp; 55 CeedInt nqpt; 56 ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr); 57 CeedInt ndof; 58 ierr = CeedBasisGetNumNodes(basis, &ndof); CeedChk(ierr); 59 const CeedInt add = (tmode == CEED_TRANSPOSE); 60 61 if (nelem != 1) 62 return CeedError(basis->ceed, 1, 63 "This backend does not support BasisApply for multiple elements"); 64 65 // Clear v if operating in transpose 66 if (tmode == CEED_TRANSPOSE) { 67 const CeedInt vsize = ncomp*ndof; 68 for (CeedInt i = 0; i < vsize; i++) 69 v[i] = (CeedScalar) 0; 70 } 71 // Tensor basis 72 if (basis->tensorbasis) { 73 switch (emode) { 74 // Interpolate to/from quadrature points 75 case CEED_EVAL_INTERP: { 76 CeedInt P = basis->P1d, Q = basis->Q1d; 77 if (tmode == CEED_TRANSPOSE) { 78 P = basis->Q1d; Q = basis->P1d; 79 } 80 CeedInt pre = ncomp*CeedIntPow(P, dim-1), post = 1; 81 CeedScalar tmp[2][ncomp*Q*CeedIntPow(P>Q?P:Q, dim-1)]; 82 for (CeedInt d=0; d<dim; d++) { 83 ierr = CeedTensorContract_Ref(basis->ceed, pre, P, post, Q, basis->interp1d, 84 tmode, add&&(d==dim-1), 85 d==0?u:tmp[d%2], d==dim-1?v:tmp[(d+1)%2]); 86 CeedChk(ierr); 87 pre /= P; 88 post *= Q; 89 } 90 } break; 91 // Evaluate the gradient to/from quadrature points 92 case CEED_EVAL_GRAD: { 93 CeedInt P = basis->P1d, Q = basis->Q1d; 94 // In CEED_NOTRANSPOSE mode: 95 // u has shape [dim, ncomp, P^dim, nelem], row-major layout 96 // v has shape [dim, ncomp, Q^dim, nelem], row-major layout 97 // In CEED_TRANSPOSE mode, the sizes of u and v are switched. 98 if (tmode == CEED_TRANSPOSE) { 99 P = basis->Q1d, Q = basis->P1d; 100 } 101 CeedScalar tmp[2][ncomp*Q*CeedIntPow(P>Q?P:Q, dim-1)]; 102 for (CeedInt p = 0; p < dim; p++) { 103 CeedInt pre = ncomp*CeedIntPow(P, dim-1), post = 1; 104 for (CeedInt d=0; d<dim; d++) { 105 ierr = CeedTensorContract_Ref(basis->ceed, pre, P, post, Q, 106 (p==d)?basis->grad1d:basis->interp1d, 107 tmode, add&&(d==dim-1), 108 (d == 0 109 ? (tmode==CEED_NOTRANSPOSE?u:u+p*ncomp*nqpt) 110 : tmp[d%2]), 111 (d == dim-1 112 ? (tmode==CEED_TRANSPOSE?v:v+p*ncomp*nqpt) 113 : tmp[(d+1)%2])); 114 CeedChk(ierr); 115 pre /= P; 116 post *= Q; 117 } 118 } 119 } break; 120 // Retrieve interpolation weights 121 case CEED_EVAL_WEIGHT: { 122 if (tmode == CEED_TRANSPOSE) 123 return CeedError(basis->ceed, 1, 124 "CEED_EVAL_WEIGHT incompatible with CEED_TRANSPOSE"); 125 CeedInt Q = basis->Q1d; 126 for (CeedInt d=0; d<dim; d++) { 127 CeedInt pre = CeedIntPow(Q, dim-d-1), post = CeedIntPow(Q, d); 128 for (CeedInt i=0; i<pre; i++) 129 for (CeedInt j=0; j<Q; j++) 130 for (CeedInt k=0; k<post; k++) 131 v[(i*Q + j)*post + k] = basis->qweight1d[j] 132 * (d == 0 ? 1 : v[(i*Q + j)*post + k]); 133 } 134 } break; 135 // Evaluate the divergence to/from the quadrature points 136 case CEED_EVAL_DIV: 137 return CeedError(basis->ceed, 1, "CEED_EVAL_DIV not supported"); 138 // Evaluate the curl to/from the quadrature points 139 case CEED_EVAL_CURL: 140 return CeedError(basis->ceed, 1, "CEED_EVAL_CURL not supported"); 141 // Take no action, BasisApply should not have been called 142 case CEED_EVAL_NONE: 143 return CeedError(basis->ceed, 1, 144 "CEED_EVAL_NONE does not make sense in this context"); 145 } 146 } else { 147 // Non-tensor basis 148 switch (emode) { 149 case CEED_EVAL_INTERP: { 150 CeedInt P = basis->P, Q = basis->Q; 151 if (tmode == CEED_TRANSPOSE) { 152 P = basis->Q; Q = basis->P; 153 } 154 ierr = CeedTensorContract_Ref(basis->ceed, ncomp, P, 1, Q, basis->interp1d, 155 tmode, add, u, v); 156 CeedChk(ierr); 157 } 158 break; 159 case CEED_EVAL_GRAD: { 160 CeedInt P = basis->P, Q = dim*basis->Q; 161 if (tmode == CEED_TRANSPOSE) { 162 P = dim*basis->Q; Q = basis->P; 163 } 164 ierr = CeedTensorContract_Ref(basis->ceed, ncomp, P, 1, Q, basis->grad1d, 165 tmode, add, u, v); 166 CeedChk(ierr); 167 } 168 break; 169 case CEED_EVAL_WEIGHT: { 170 if (tmode == CEED_TRANSPOSE) 171 return CeedError(basis->ceed, 1, 172 "CEED_EVAL_WEIGHT incompatible with CEED_TRANSPOSE"); 173 for (CeedInt i=0; i<nqpt; i++) 174 v[i] = basis->qweight1d[i]; 175 } break; 176 case CEED_EVAL_DIV: 177 return CeedError(basis->ceed, 1, "CEED_EVAL_DIV not supported"); 178 case CEED_EVAL_CURL: 179 return CeedError(basis->ceed, 1, "CEED_EVAL_CURL not supported"); 180 case CEED_EVAL_NONE: 181 return CeedError(basis->ceed, 1, 182 "CEED_EVAL_NONE does not make sense in this context"); 183 } 184 } 185 return 0; 186 } 187 188 static int CeedBasisDestroy_Ref(CeedBasis basis) { 189 return 0; 190 } 191 192 int CeedBasisCreateTensorH1_Ref(Ceed ceed, CeedInt dim, CeedInt P1d, 193 CeedInt Q1d, const CeedScalar *interp1d, 194 const CeedScalar *grad1d, 195 const CeedScalar *qref1d, 196 const CeedScalar *qweight1d, 197 CeedBasis basis) { 198 basis->Apply = CeedBasisApply_Ref; 199 basis->Destroy = CeedBasisDestroy_Ref; 200 return 0; 201 } 202 203 int CeedBasisCreateH1_Ref(Ceed ceed, CeedElemTopology topo, CeedInt dim, 204 CeedInt ndof, CeedInt nqpts, 205 const CeedScalar *interp, 206 const CeedScalar *grad, 207 const CeedScalar *qref, 208 const CeedScalar *qweight, 209 CeedBasis basis) { 210 basis->Apply = CeedBasisApply_Ref; 211 basis->Destroy = CeedBasisDestroy_Ref; 212 return 0; 213 } 214