1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ed264d09SValeria Barra // 4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5*3d8e8822SJeremy L Thompson // 6*3d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed element and spectral 7ed264d09SValeria Barra // element discretizations for exascale applications. For more information and 8*3d8e8822SJeremy L Thompson // source code availability see http://github.com/ceed 9ed264d09SValeria Barra // 10ed264d09SValeria Barra // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11ed264d09SValeria Barra // a collaborative effort of two U.S. Department of Energy organizations (Office 12ed264d09SValeria Barra // of Science and the National Nuclear Security Administration) responsible for 13ed264d09SValeria Barra // the planning and preparation of a capable exascale ecosystem, including 14ed264d09SValeria Barra // software, applications, hardware, advanced system engineering and early 15ed264d09SValeria Barra // testbed platforms, in support of the nation's exascale computing imperative. 16ed264d09SValeria Barra 17ed264d09SValeria Barra /// @file 18ed264d09SValeria Barra /// libCEED QFunctions for diffusion operator example for a scalar field on the sphere using PETSc 19ed264d09SValeria Barra 20f6b55d2cSvaleriabarra #ifndef bp3sphere_h 21f6b55d2cSvaleriabarra #define bp3sphere_h 22f6b55d2cSvaleriabarra 23ed264d09SValeria Barra #include <math.h> 24ed264d09SValeria Barra 25e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 26ed264d09SValeria Barra // This QFunction sets up the geometric factors required for integration and 27ed264d09SValeria Barra // coordinate transformations when reference coordinates have a different 28ed264d09SValeria Barra // dimension than the one of physical coordinates 29ed264d09SValeria Barra // 30ed264d09SValeria Barra // Reference (parent) 2D coordinates: X \in [-1, 1]^2 31ed264d09SValeria Barra // 32ed264d09SValeria Barra // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 33ed264d09SValeria Barra // with R radius of the sphere 34ed264d09SValeria Barra // 35ed264d09SValeria Barra // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3 36ed264d09SValeria Barra // with l half edge of the cube inscribed in the sphere 37ed264d09SValeria Barra // 38ed264d09SValeria Barra // Change of coordinates matrix computed by the library: 39ed264d09SValeria Barra // (physical 3D coords relative to reference 2D coords) 40ed264d09SValeria Barra // dxx_j/dX_i (indicial notation) [3 * 2] 41ed264d09SValeria Barra // 42ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D): 43ed264d09SValeria Barra // dx_i/dxx_j (indicial notation) [3 * 3] 44ed264d09SValeria Barra // 45ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to X (reference 2D): 46ed264d09SValeria Barra // (by chain rule) 47ed264d09SValeria Barra // dx_i/dX_j [3 * 2] = dx_i/dxx_k [3 * 3] * dxx_k/dX_j [3 * 2] 48ed264d09SValeria Barra // 499b072555Sjeremylt // mod_J is given by the magnitude of the cross product of the columns of dx_i/dX_j 50ed264d09SValeria Barra // 519b072555Sjeremylt // The quadrature data is stored in the array q_data. 52ed264d09SValeria Barra // 53ed264d09SValeria Barra // We require the determinant of the Jacobian to properly compute integrals of 54ed264d09SValeria Barra // the form: int( u v ) 55ed264d09SValeria Barra // 569b072555Sjeremylt // q_data[0]: mod_J * w 57ed264d09SValeria Barra // 58ed264d09SValeria Barra // We use the Moore–Penrose (left) pseudoinverse of dx_i/dX_j, to compute dX_i/dx_j (and its transpose), 59ed264d09SValeria Barra // needed to properly compute integrals of the form: int( gradv gradu ) 60ed264d09SValeria Barra // 61ed264d09SValeria Barra // dX_i/dx_j [2 * 3] = (dx_i/dX_j)+ = (dxdX^T dxdX)^(-1) dxdX 62ed264d09SValeria Barra // 63ac4340cfSJed Brown // and the product simplifies to yield the contravariant metric tensor 64ac4340cfSJed Brown // 65ac4340cfSJed Brown // g^{ij} = dX_i/dx_k dX_j/dx_k = (dxdX^T dxdX)^{-1} 66ac4340cfSJed Brown // 6708fade8cSvaleriabarra // Stored: g^{ij} (in Voigt convention) in 6808fade8cSvaleriabarra // 699b072555Sjeremylt // q_data[1:3]: [dXdxdXdxT00 dXdxdXdxT01] 7008fade8cSvaleriabarra // [dXdxdXdxT01 dXdxdXdxT11] 71ed264d09SValeria Barra // ----------------------------------------------------------------------------- 72ed264d09SValeria Barra CEED_QFUNCTION(SetupDiffGeo)(void *ctx, CeedInt Q, 73ed264d09SValeria Barra const CeedScalar *const *in, 74ed264d09SValeria Barra CeedScalar *const *out) { 75ed264d09SValeria Barra const CeedScalar *X = in[0], *J = in[1], *w = in[2]; 769b072555Sjeremylt CeedScalar *q_data = out[0]; 77ed264d09SValeria Barra 78ed264d09SValeria Barra // Quadrature Point Loop 79ed264d09SValeria Barra CeedPragmaSIMD 80ed264d09SValeria Barra for (CeedInt i=0; i<Q; i++) { 81ed264d09SValeria Barra // Read global Cartesian coordinates 82ed264d09SValeria Barra const CeedScalar xx[3] = {X[i+0*Q], 83ed264d09SValeria Barra X[i+1*Q], 84ed264d09SValeria Barra X[i+2*Q] 85ed264d09SValeria Barra }; 86ed264d09SValeria Barra 87ed264d09SValeria Barra // Read dxxdX Jacobian entries, stored as 88ed264d09SValeria Barra // 0 3 89ed264d09SValeria Barra // 1 4 90ed264d09SValeria Barra // 2 5 91ed264d09SValeria Barra const CeedScalar dxxdX[3][2] = {{J[i+Q*0], 92ed264d09SValeria Barra J[i+Q*3]}, 93ed264d09SValeria Barra {J[i+Q*1], 94ed264d09SValeria Barra J[i+Q*4]}, 95ed264d09SValeria Barra {J[i+Q*2], 96ed264d09SValeria Barra J[i+Q*5]} 97ed264d09SValeria Barra }; 98ed264d09SValeria Barra 99ed264d09SValeria Barra // Setup 100ed264d09SValeria Barra // x = xx (xx^T xx)^{-1/2} 101ed264d09SValeria Barra // dx/dxx = I (xx^T xx)^{-1/2} - xx xx^T (xx^T xx)^{-3/2} 1029b072555Sjeremylt const CeedScalar mod_xx_sq = xx[0]*xx[0]+xx[1]*xx[1]+xx[2]*xx[2]; 1039b072555Sjeremylt CeedScalar xx_sq[3][3]; 104ed264d09SValeria Barra for (int j=0; j<3; j++) 105ed264d09SValeria Barra for (int k=0; k<3; k++) 1069b072555Sjeremylt xx_sq[j][k] = xx[j]*xx[k] / (sqrt(mod_xx_sq) * mod_xx_sq); 107ed264d09SValeria Barra 1089b072555Sjeremylt const CeedScalar dxdxx[3][3] = {{1./sqrt(mod_xx_sq) - xx_sq[0][0], 1099b072555Sjeremylt -xx_sq[0][1], 1109b072555Sjeremylt -xx_sq[0][2]}, 1119b072555Sjeremylt {-xx_sq[1][0], 1129b072555Sjeremylt 1./sqrt(mod_xx_sq) - xx_sq[1][1], 1139b072555Sjeremylt -xx_sq[1][2]}, 1149b072555Sjeremylt {-xx_sq[2][0], 1159b072555Sjeremylt -xx_sq[2][1], 1169b072555Sjeremylt 1./sqrt(mod_xx_sq) - xx_sq[2][2]} 117ed264d09SValeria Barra }; 118ed264d09SValeria Barra 119ed264d09SValeria Barra CeedScalar dxdX[3][2]; 120ed264d09SValeria Barra for (int j=0; j<3; j++) 121ed264d09SValeria Barra for (int k=0; k<2; k++) { 122ed264d09SValeria Barra dxdX[j][k] = 0; 123ed264d09SValeria Barra for (int l=0; l<3; l++) 124ed264d09SValeria Barra dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k]; 125ed264d09SValeria Barra } 126ed264d09SValeria Barra 127ed264d09SValeria Barra // J is given by the cross product of the columns of dxdX 128ed264d09SValeria Barra const CeedScalar J[3]= {dxdX[1][0]*dxdX[2][1] - dxdX[2][0]*dxdX[1][1], 129ed264d09SValeria Barra dxdX[2][0]*dxdX[0][1] - dxdX[0][0]*dxdX[2][1], 130ed264d09SValeria Barra dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1] 131ed264d09SValeria Barra }; 132ed264d09SValeria Barra 133ed264d09SValeria Barra // Use the magnitude of J as our detJ (volume scaling factor) 1349b072555Sjeremylt const CeedScalar mod_J = sqrt(J[0]*J[0]+J[1]*J[1]+J[2]*J[2]); 135ed264d09SValeria Barra 1369b072555Sjeremylt // Interp-to-Interp q_data 1379b072555Sjeremylt q_data[i+Q*0] = mod_J * w[i]; 138ed264d09SValeria Barra 13908fade8cSvaleriabarra // dxdX_k,j * dxdX_j,k 140ed264d09SValeria Barra CeedScalar dxdXTdxdX[2][2]; 141ed264d09SValeria Barra for (int j=0; j<2; j++) 142ed264d09SValeria Barra for (int k=0; k<2; k++) { 143ed264d09SValeria Barra dxdXTdxdX[j][k] = 0; 144ed264d09SValeria Barra for (int l=0; l<3; l++) 145ed264d09SValeria Barra dxdXTdxdX[j][k] += dxdX[l][j]*dxdX[l][k]; 146ed264d09SValeria Barra } 147ed264d09SValeria Barra 148ed264d09SValeria Barra const CeedScalar detdxdXTdxdX = dxdXTdxdX[0][0] * dxdXTdxdX[1][1] 149ed264d09SValeria Barra -dxdXTdxdX[1][0] * dxdXTdxdX[0][1]; 150ed264d09SValeria Barra 15108fade8cSvaleriabarra // Compute inverse of dxdXTdxdX, which is the 2x2 contravariant metric tensor g^{ij} 1529b072555Sjeremylt CeedScalar dxdXTdxdX_inv[2][2]; 1539b072555Sjeremylt dxdXTdxdX_inv[0][0] = dxdXTdxdX[1][1] / detdxdXTdxdX; 1549b072555Sjeremylt dxdXTdxdX_inv[0][1] = -dxdXTdxdX[0][1] / detdxdXTdxdX; 1559b072555Sjeremylt dxdXTdxdX_inv[1][0] = -dxdXTdxdX[1][0] / detdxdXTdxdX; 1569b072555Sjeremylt dxdXTdxdX_inv[1][1] = dxdXTdxdX[0][0] / detdxdXTdxdX; 157ed264d09SValeria Barra 158ed264d09SValeria Barra // Stored in Voigt convention 1599b072555Sjeremylt q_data[i+Q*1] = dxdXTdxdX_inv[0][0]; 1609b072555Sjeremylt q_data[i+Q*2] = dxdXTdxdX_inv[1][1]; 1619b072555Sjeremylt q_data[i+Q*3] = dxdXTdxdX_inv[0][1]; 162ed264d09SValeria Barra } // End of Quadrature Point Loop 163ed264d09SValeria Barra 164ed264d09SValeria Barra // Return 165ed264d09SValeria Barra return 0; 166ed264d09SValeria Barra } 167ed264d09SValeria Barra 168e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 169ed264d09SValeria Barra // This QFunction sets up the rhs and true solution for the problem 170ed264d09SValeria Barra // ----------------------------------------------------------------------------- 171ed264d09SValeria Barra CEED_QFUNCTION(SetupDiffRhs)(void *ctx, CeedInt Q, 172ed264d09SValeria Barra const CeedScalar *const *in, 173ed264d09SValeria Barra CeedScalar *const *out) { 174ed264d09SValeria Barra // Inputs 1759b072555Sjeremylt const CeedScalar *X = in[0], *q_data = in[1]; 176ed264d09SValeria Barra // Outputs 177ed264d09SValeria Barra CeedScalar *true_soln = out[0], *rhs = out[1]; 178ed264d09SValeria Barra 179ed264d09SValeria Barra // Context 180ed264d09SValeria Barra const CeedScalar *context = (const CeedScalar*)ctx; 181ed264d09SValeria Barra const CeedScalar R = context[0]; 182ed264d09SValeria Barra 183ed264d09SValeria Barra // Quadrature Point Loop 184ed264d09SValeria Barra CeedPragmaSIMD 185ed264d09SValeria Barra for (CeedInt i=0; i<Q; i++) { 186ed264d09SValeria Barra // Read global Cartesian coordinates 187ed264d09SValeria Barra CeedScalar x = X[i+Q*0], y = X[i+Q*1], z = X[i+Q*2]; 188ed264d09SValeria Barra // Normalize quadrature point coordinates to sphere 189ed264d09SValeria Barra CeedScalar rad = sqrt(x*x + y*y + z*z); 190ed264d09SValeria Barra x *= R / rad; 191ed264d09SValeria Barra y *= R / rad; 192ed264d09SValeria Barra z *= R / rad; 193ed264d09SValeria Barra // Compute latitude and longitude 194ed264d09SValeria Barra const CeedScalar theta = asin(z / R); // latitude 195ed264d09SValeria Barra const CeedScalar lambda = atan2(y, x); // longitude 196ed264d09SValeria Barra 197ed264d09SValeria Barra true_soln[i+Q*0] = sin(lambda) * cos(theta); 198ed264d09SValeria Barra 1999b072555Sjeremylt rhs[i+Q*0] = q_data[i+Q*0] * 2 * sin(lambda)*cos(theta) / (R*R); 200ed264d09SValeria Barra 201ed264d09SValeria Barra } // End of Quadrature Point Loop 202ed264d09SValeria Barra 203ed264d09SValeria Barra return 0; 204ed264d09SValeria Barra } 205ed264d09SValeria Barra 206e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 207ed264d09SValeria Barra // This QFunction applies the diffusion operator for a scalar field. 208ed264d09SValeria Barra // 209ed264d09SValeria Barra // Inputs: 210ed264d09SValeria Barra // ug - Input vector gradient at quadrature points 2119b072555Sjeremylt // q_data - Geometric factors 212ed264d09SValeria Barra // 213ed264d09SValeria Barra // Output: 214ed264d09SValeria Barra // vg - Output vector (test functions) gradient at quadrature points 215ed264d09SValeria Barra // 216ed264d09SValeria Barra // ----------------------------------------------------------------------------- 217ed264d09SValeria Barra CEED_QFUNCTION(Diff)(void *ctx, CeedInt Q, 218ed264d09SValeria Barra const CeedScalar *const *in, CeedScalar *const *out) { 219ed264d09SValeria Barra // Inputs 2209b072555Sjeremylt const CeedScalar *ug = in[0], *q_data = in[1]; 221ed264d09SValeria Barra // Outputs 222ed264d09SValeria Barra CeedScalar *vg = out[0]; 223ed264d09SValeria Barra 224ed264d09SValeria Barra // Quadrature Point Loop 225ed264d09SValeria Barra CeedPragmaSIMD 226ed264d09SValeria Barra for (CeedInt i=0; i<Q; i++) { 227ed264d09SValeria Barra // Read spatial derivatives of u 228ed264d09SValeria Barra const CeedScalar du[2] = {ug[i+Q*0], 229ed264d09SValeria Barra ug[i+Q*1] 230ed264d09SValeria Barra }; 2319b072555Sjeremylt // Read q_data 2329b072555Sjeremylt const CeedScalar w_det_J = q_data[i+Q*0]; 2339b072555Sjeremylt // -- Grad-to-Grad q_data 234ed264d09SValeria Barra // ---- dXdx_j,k * dXdx_k,j 2359b072555Sjeremylt const CeedScalar dXdxdXdx_T[2][2] = {{q_data[i+Q*1], 2369b072555Sjeremylt q_data[i+Q*3]}, 2379b072555Sjeremylt {q_data[i+Q*3], 2389b072555Sjeremylt q_data[i+Q*2]} 239ed264d09SValeria Barra }; 240ed264d09SValeria Barra 241ed264d09SValeria Barra for (int j=0; j<2; j++) // j = direction of vg 2429b072555Sjeremylt vg[i+j*Q] = w_det_J * (du[0] * dXdxdXdx_T[0][j] + 2439b072555Sjeremylt du[1] * dXdxdXdx_T[1][j]); 244ed264d09SValeria Barra 245ed264d09SValeria Barra } // End of Quadrature Point Loop 246ed264d09SValeria Barra 247ed264d09SValeria Barra return 0; 248ed264d09SValeria Barra } 249ed264d09SValeria Barra // ----------------------------------------------------------------------------- 250f6b55d2cSvaleriabarra 251f6b55d2cSvaleriabarra #endif // bp3sphere_h 252