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