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