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 /// libCEED QFunctions for diffusion operator example for a scalar field on the sphere using PETSc 19 20 #ifndef bp3sphere_h 21 #define bp3sphere_h 22 #include <ceed.h> 23 24 #ifndef __CUDACC__ 25 # include <math.h> 26 #endif 27 28 // ***************************************************************************** 29 // This QFunction sets up the geometric factors required for integration and 30 // coordinate transformations when reference coordinates have a different 31 // dimension than the one of physical coordinates 32 // 33 // Reference (parent) 2D coordinates: X \in [-1, 1]^2 34 // 35 // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 36 // with R radius of the sphere 37 // 38 // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3 39 // with l half edge of the cube inscribed in the sphere 40 // 41 // Change of coordinates matrix computed by the library: 42 // (physical 3D coords relative to reference 2D coords) 43 // dxx_j/dX_i (indicial notation) [3 * 2] 44 // 45 // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D): 46 // dx_i/dxx_j (indicial notation) [3 * 3] 47 // 48 // Change of coordinates x (on the 2D manifold) relative to X (reference 2D): 49 // (by chain rule) 50 // dx_i/dX_j [3 * 2] = dx_i/dxx_k [3 * 3] * dxx_k/dX_j [3 * 2] 51 // 52 // modJ is given by the magnitude of the cross product of the columns of dx_i/dX_j 53 // 54 // The quadrature data is stored in the array qdata. 55 // 56 // We require the determinant of the Jacobian to properly compute integrals of 57 // the form: int( u v ) 58 // 59 // qdata[0]: modJ * w 60 // 61 // We use the Moore–Penrose (left) pseudoinverse of dx_i/dX_j, to compute dX_i/dx_j (and its transpose), 62 // needed to properly compute integrals of the form: int( gradv gradu ) 63 // 64 // dX_i/dx_j [2 * 3] = (dx_i/dX_j)+ = (dxdX^T dxdX)^(-1) dxdX 65 // 66 // and the product simplifies to yield the contravariant metric tensor 67 // 68 // g^{ij} = dX_i/dx_k dX_j/dx_k = (dxdX^T dxdX)^{-1} 69 // 70 // Stored: g^{ij} (in Voigt convention) in 71 // 72 // qdata[1:3]: [dXdxdXdxT00 dXdxdXdxT01] 73 // [dXdxdXdxT01 dXdxdXdxT11] 74 // ***************************************************************************** 75 76 // ----------------------------------------------------------------------------- 77 CEED_QFUNCTION(SetupDiffGeo)(void *ctx, CeedInt Q, 78 const CeedScalar *const *in, 79 CeedScalar *const *out) { 80 const CeedScalar *X = in[0], *J = in[1], *w = in[2]; 81 CeedScalar *qdata = out[0]; 82 83 // Quadrature Point Loop 84 CeedPragmaSIMD 85 for (CeedInt i=0; i<Q; i++) { 86 // Read global Cartesian coordinates 87 const CeedScalar xx[3] = {X[i+0*Q], 88 X[i+1*Q], 89 X[i+2*Q] 90 }; 91 92 // Read dxxdX Jacobian entries, stored as 93 // 0 3 94 // 1 4 95 // 2 5 96 const CeedScalar dxxdX[3][2] = {{J[i+Q*0], 97 J[i+Q*3]}, 98 {J[i+Q*1], 99 J[i+Q*4]}, 100 {J[i+Q*2], 101 J[i+Q*5]} 102 }; 103 104 // Setup 105 // x = xx (xx^T xx)^{-1/2} 106 // dx/dxx = I (xx^T xx)^{-1/2} - xx xx^T (xx^T xx)^{-3/2} 107 const CeedScalar modxxsq = xx[0]*xx[0]+xx[1]*xx[1]+xx[2]*xx[2]; 108 CeedScalar xxsq[3][3]; 109 for (int j=0; j<3; j++) 110 for (int k=0; k<3; k++) 111 xxsq[j][k] = xx[j]*xx[k] / (sqrt(modxxsq) * modxxsq); 112 113 const CeedScalar dxdxx[3][3] = {{1./sqrt(modxxsq) - xxsq[0][0], 114 -xxsq[0][1], 115 -xxsq[0][2]}, 116 {-xxsq[1][0], 117 1./sqrt(modxxsq) - xxsq[1][1], 118 -xxsq[1][2]}, 119 {-xxsq[2][0], 120 -xxsq[2][1], 121 1./sqrt(modxxsq) - xxsq[2][2]} 122 }; 123 124 CeedScalar dxdX[3][2]; 125 for (int j=0; j<3; j++) 126 for (int k=0; k<2; k++) { 127 dxdX[j][k] = 0; 128 for (int l=0; l<3; l++) 129 dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k]; 130 } 131 132 // J is given by the cross product of the columns of dxdX 133 const CeedScalar J[3]= {dxdX[1][0]*dxdX[2][1] - dxdX[2][0]*dxdX[1][1], 134 dxdX[2][0]*dxdX[0][1] - dxdX[0][0]*dxdX[2][1], 135 dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1] 136 }; 137 138 // Use the magnitude of J as our detJ (volume scaling factor) 139 const CeedScalar modJ = sqrt(J[0]*J[0]+J[1]*J[1]+J[2]*J[2]); 140 141 // Interp-to-Interp qdata 142 qdata[i+Q*0] = modJ * w[i]; 143 144 // dxdX_k,j * dxdX_j,k 145 CeedScalar dxdXTdxdX[2][2]; 146 for (int j=0; j<2; j++) 147 for (int k=0; k<2; k++) { 148 dxdXTdxdX[j][k] = 0; 149 for (int l=0; l<3; l++) 150 dxdXTdxdX[j][k] += dxdX[l][j]*dxdX[l][k]; 151 } 152 153 const CeedScalar detdxdXTdxdX = dxdXTdxdX[0][0] * dxdXTdxdX[1][1] 154 -dxdXTdxdX[1][0] * dxdXTdxdX[0][1]; 155 156 // Compute inverse of dxdXTdxdX, which is the 2x2 contravariant metric tensor g^{ij} 157 CeedScalar dxdXTdxdXinv[2][2]; 158 dxdXTdxdXinv[0][0] = dxdXTdxdX[1][1] / detdxdXTdxdX; 159 dxdXTdxdXinv[0][1] = -dxdXTdxdX[0][1] / detdxdXTdxdX; 160 dxdXTdxdXinv[1][0] = -dxdXTdxdX[1][0] / detdxdXTdxdX; 161 dxdXTdxdXinv[1][1] = dxdXTdxdX[0][0] / detdxdXTdxdX; 162 163 // Stored in Voigt convention 164 qdata[i+Q*1] = dxdXTdxdXinv[0][0]; 165 qdata[i+Q*2] = dxdXTdxdXinv[1][1]; 166 qdata[i+Q*3] = dxdXTdxdXinv[0][1]; 167 } // End of Quadrature Point Loop 168 169 // Return 170 return 0; 171 } 172 173 // ***************************************************************************** 174 // This QFunction sets up the rhs and true solution for the problem 175 // ***************************************************************************** 176 177 // ----------------------------------------------------------------------------- 178 CEED_QFUNCTION(SetupDiffRhs)(void *ctx, CeedInt Q, 179 const CeedScalar *const *in, 180 CeedScalar *const *out) { 181 // Inputs 182 const CeedScalar *X = in[0], *qdata = in[1]; 183 // Outputs 184 CeedScalar *true_soln = out[0], *rhs = out[1]; 185 186 // Context 187 const CeedScalar *context = (const CeedScalar*)ctx; 188 const CeedScalar R = context[0]; 189 190 // Quadrature Point Loop 191 CeedPragmaSIMD 192 for (CeedInt i=0; i<Q; i++) { 193 // Read global Cartesian coordinates 194 CeedScalar x = X[i+Q*0], y = X[i+Q*1], z = X[i+Q*2]; 195 // Normalize quadrature point coordinates to sphere 196 CeedScalar rad = sqrt(x*x + y*y + z*z); 197 x *= R / rad; 198 y *= R / rad; 199 z *= R / rad; 200 // Compute latitude and longitude 201 const CeedScalar theta = asin(z / R); // latitude 202 const CeedScalar lambda = atan2(y, x); // longitude 203 204 true_soln[i+Q*0] = sin(lambda) * cos(theta); 205 206 rhs[i+Q*0] = qdata[i+Q*0] * 2 * sin(lambda)*cos(theta) / (R*R); 207 208 } // End of Quadrature Point Loop 209 210 return 0; 211 } 212 213 // ***************************************************************************** 214 // This QFunction applies the diffusion operator for a scalar field. 215 // 216 // Inputs: 217 // ug - Input vector gradient at quadrature points 218 // qdata - Geometric factors 219 // 220 // Output: 221 // vg - Output vector (test functions) gradient at quadrature points 222 // 223 // ***************************************************************************** 224 225 // ----------------------------------------------------------------------------- 226 CEED_QFUNCTION(Diff)(void *ctx, CeedInt Q, 227 const CeedScalar *const *in, CeedScalar *const *out) { 228 // Inputs 229 const CeedScalar *ug = in[0], *qdata = in[1]; 230 // Outputs 231 CeedScalar *vg = out[0]; 232 233 // Quadrature Point Loop 234 CeedPragmaSIMD 235 for (CeedInt i=0; i<Q; i++) { 236 // Read spatial derivatives of u 237 const CeedScalar du[2] = {ug[i+Q*0], 238 ug[i+Q*1] 239 }; 240 // Read qdata 241 const CeedScalar wJ = qdata[i+Q*0]; 242 // -- Grad-to-Grad qdata 243 // ---- dXdx_j,k * dXdx_k,j 244 const CeedScalar dXdxdXdxT[2][2] = {{qdata[i+Q*1], 245 qdata[i+Q*3]}, 246 {qdata[i+Q*3], 247 qdata[i+Q*2]} 248 }; 249 250 for (int j=0; j<2; j++) // j = direction of vg 251 vg[i+j*Q] = wJ * (du[0] * dXdxdXdxT[0][j] + 252 du[1] * dXdxdXdxT[1][j]); 253 254 } // End of Quadrature Point Loop 255 256 return 0; 257 } 258 // ----------------------------------------------------------------------------- 259 260 #endif // bp3sphere_h 261