13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ed264d09SValeria Barra // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 53d8e8822SJeremy L Thompson // 63d8e8822SJeremy 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 83d8e8822SJeremy 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 23c9c2c079SJeremy L Thompson #include <ceed.h> 24ed264d09SValeria Barra #include <math.h> 25ed264d09SValeria Barra 26e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 27ed264d09SValeria Barra // This QFunction sets up the geometric factors required for integration and 28ed264d09SValeria Barra // coordinate transformations when reference coordinates have a different 29ed264d09SValeria Barra // dimension than the one of physical coordinates 30ed264d09SValeria Barra // 31ed264d09SValeria Barra // Reference (parent) 2D coordinates: X \in [-1, 1]^2 32ed264d09SValeria Barra // 33ed264d09SValeria Barra // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 34ed264d09SValeria Barra // with R radius of the sphere 35ed264d09SValeria Barra // 36ed264d09SValeria Barra // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3 37ed264d09SValeria Barra // with l half edge of the cube inscribed in the sphere 38ed264d09SValeria Barra // 39ed264d09SValeria Barra // Change of coordinates matrix computed by the library: 40ed264d09SValeria Barra // (physical 3D coords relative to reference 2D coords) 41ed264d09SValeria Barra // dxx_j/dX_i (indicial notation) [3 * 2] 42ed264d09SValeria Barra // 43ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D): 44ed264d09SValeria Barra // dx_i/dxx_j (indicial notation) [3 * 3] 45ed264d09SValeria Barra // 46ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to X (reference 2D): 47ed264d09SValeria Barra // (by chain rule) 48ed264d09SValeria Barra // dx_i/dX_j [3 * 2] = dx_i/dxx_k [3 * 3] * dxx_k/dX_j [3 * 2] 49ed264d09SValeria Barra // 509b072555Sjeremylt // mod_J is given by the magnitude of the cross product of the columns of dx_i/dX_j 51ed264d09SValeria Barra // 529b072555Sjeremylt // The quadrature data is stored in the array q_data. 53ed264d09SValeria Barra // 54ed264d09SValeria Barra // We require the determinant of the Jacobian to properly compute integrals of 55ed264d09SValeria Barra // the form: int( u v ) 56ed264d09SValeria Barra // 579b072555Sjeremylt // q_data[0]: mod_J * w 58ed264d09SValeria Barra // 59ed264d09SValeria Barra // We use the Moore–Penrose (left) pseudoinverse of dx_i/dX_j, to compute dX_i/dx_j (and its transpose), 60ed264d09SValeria Barra // needed to properly compute integrals of the form: int( gradv gradu ) 61ed264d09SValeria Barra // 62ed264d09SValeria Barra // dX_i/dx_j [2 * 3] = (dx_i/dX_j)+ = (dxdX^T dxdX)^(-1) dxdX 63ed264d09SValeria Barra // 64ac4340cfSJed Brown // and the product simplifies to yield the contravariant metric tensor 65ac4340cfSJed Brown // 66ac4340cfSJed Brown // g^{ij} = dX_i/dx_k dX_j/dx_k = (dxdX^T dxdX)^{-1} 67ac4340cfSJed Brown // 6808fade8cSvaleriabarra // Stored: g^{ij} (in Voigt convention) in 6908fade8cSvaleriabarra // 709b072555Sjeremylt // q_data[1:3]: [dXdxdXdxT00 dXdxdXdxT01] 7108fade8cSvaleriabarra // [dXdxdXdxT01 dXdxdXdxT11] 72ed264d09SValeria Barra // ----------------------------------------------------------------------------- 73*2b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffGeo)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 74ed264d09SValeria Barra const CeedScalar *X = in[0], *J = in[1], *w = in[2]; 759b072555Sjeremylt CeedScalar *q_data = out[0]; 76ed264d09SValeria Barra 77ed264d09SValeria Barra // Quadrature Point Loop 78*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 79ed264d09SValeria Barra // Read global Cartesian coordinates 80*2b730f8bSJeremy L Thompson const CeedScalar xx[3] = {X[i + 0 * Q], X[i + 1 * Q], X[i + 2 * Q]}; 81ed264d09SValeria Barra 82ed264d09SValeria Barra // Read dxxdX Jacobian entries, stored as 83ed264d09SValeria Barra // 0 3 84ed264d09SValeria Barra // 1 4 85ed264d09SValeria Barra // 2 5 86*2b730f8bSJeremy L Thompson const CeedScalar dxxdX[3][2] = { 87*2b730f8bSJeremy L Thompson {J[i + Q * 0], J[i + Q * 3]}, 88*2b730f8bSJeremy L Thompson {J[i + Q * 1], J[i + Q * 4]}, 89*2b730f8bSJeremy L Thompson {J[i + Q * 2], J[i + Q * 5]} 90ed264d09SValeria Barra }; 91ed264d09SValeria Barra 92ed264d09SValeria Barra // Setup 93ed264d09SValeria Barra // x = xx (xx^T xx)^{-1/2} 94ed264d09SValeria Barra // dx/dxx = I (xx^T xx)^{-1/2} - xx xx^T (xx^T xx)^{-3/2} 959b072555Sjeremylt const CeedScalar mod_xx_sq = xx[0] * xx[0] + xx[1] * xx[1] + xx[2] * xx[2]; 969b072555Sjeremylt CeedScalar xx_sq[3][3]; 97*2b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) { 98*2b730f8bSJeremy L Thompson for (int k = 0; k < 3; k++) xx_sq[j][k] = xx[j] * xx[k] / (sqrt(mod_xx_sq) * mod_xx_sq); 99*2b730f8bSJeremy L Thompson } 100ed264d09SValeria Barra 101*2b730f8bSJeremy L Thompson const CeedScalar dxdxx[3][3] = { 102*2b730f8bSJeremy L Thompson {1. / sqrt(mod_xx_sq) - xx_sq[0][0], -xx_sq[0][1], -xx_sq[0][2] }, 103*2b730f8bSJeremy L Thompson {-xx_sq[1][0], 1. / sqrt(mod_xx_sq) - xx_sq[1][1], -xx_sq[1][2] }, 104*2b730f8bSJeremy L Thompson {-xx_sq[2][0], -xx_sq[2][1], 1. / sqrt(mod_xx_sq) - xx_sq[2][2]} 105ed264d09SValeria Barra }; 106ed264d09SValeria Barra 107ed264d09SValeria Barra CeedScalar dxdX[3][2]; 108*2b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) { 109ed264d09SValeria Barra for (int k = 0; k < 2; k++) { 110ed264d09SValeria Barra dxdX[j][k] = 0; 111*2b730f8bSJeremy L Thompson for (int l = 0; l < 3; l++) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k]; 112*2b730f8bSJeremy L Thompson } 113ed264d09SValeria Barra } 114ed264d09SValeria Barra 115ed264d09SValeria Barra // J is given by the cross product of the columns of dxdX 116*2b730f8bSJeremy L Thompson const CeedScalar J[3] = {dxdX[1][0] * dxdX[2][1] - dxdX[2][0] * dxdX[1][1], dxdX[2][0] * dxdX[0][1] - dxdX[0][0] * dxdX[2][1], 117*2b730f8bSJeremy L Thompson dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]}; 118ed264d09SValeria Barra 119ed264d09SValeria Barra // Use the magnitude of J as our detJ (volume scaling factor) 1209b072555Sjeremylt const CeedScalar mod_J = sqrt(J[0] * J[0] + J[1] * J[1] + J[2] * J[2]); 121ed264d09SValeria Barra 1229b072555Sjeremylt // Interp-to-Interp q_data 1239b072555Sjeremylt q_data[i + Q * 0] = mod_J * w[i]; 124ed264d09SValeria Barra 12508fade8cSvaleriabarra // dxdX_k,j * dxdX_j,k 126ed264d09SValeria Barra CeedScalar dxdXTdxdX[2][2]; 127*2b730f8bSJeremy L Thompson for (int j = 0; j < 2; j++) { 128ed264d09SValeria Barra for (int k = 0; k < 2; k++) { 129ed264d09SValeria Barra dxdXTdxdX[j][k] = 0; 130*2b730f8bSJeremy L Thompson for (int l = 0; l < 3; l++) dxdXTdxdX[j][k] += dxdX[l][j] * dxdX[l][k]; 131*2b730f8bSJeremy L Thompson } 132ed264d09SValeria Barra } 133ed264d09SValeria Barra 134*2b730f8bSJeremy L Thompson const CeedScalar detdxdXTdxdX = dxdXTdxdX[0][0] * dxdXTdxdX[1][1] - dxdXTdxdX[1][0] * dxdXTdxdX[0][1]; 135ed264d09SValeria Barra 13608fade8cSvaleriabarra // Compute inverse of dxdXTdxdX, which is the 2x2 contravariant metric tensor g^{ij} 1379b072555Sjeremylt CeedScalar dxdXTdxdX_inv[2][2]; 1389b072555Sjeremylt dxdXTdxdX_inv[0][0] = dxdXTdxdX[1][1] / detdxdXTdxdX; 1399b072555Sjeremylt dxdXTdxdX_inv[0][1] = -dxdXTdxdX[0][1] / detdxdXTdxdX; 1409b072555Sjeremylt dxdXTdxdX_inv[1][0] = -dxdXTdxdX[1][0] / detdxdXTdxdX; 1419b072555Sjeremylt dxdXTdxdX_inv[1][1] = dxdXTdxdX[0][0] / detdxdXTdxdX; 142ed264d09SValeria Barra 143ed264d09SValeria Barra // Stored in Voigt convention 1449b072555Sjeremylt q_data[i + Q * 1] = dxdXTdxdX_inv[0][0]; 1459b072555Sjeremylt q_data[i + Q * 2] = dxdXTdxdX_inv[1][1]; 1469b072555Sjeremylt q_data[i + Q * 3] = dxdXTdxdX_inv[0][1]; 147ed264d09SValeria Barra } // End of Quadrature Point Loop 148ed264d09SValeria Barra 149ed264d09SValeria Barra // Return 150ed264d09SValeria Barra return 0; 151ed264d09SValeria Barra } 152ed264d09SValeria Barra 153e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 154ed264d09SValeria Barra // This QFunction sets up the rhs and true solution for the problem 155ed264d09SValeria Barra // ----------------------------------------------------------------------------- 156*2b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffRhs)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 157ed264d09SValeria Barra // Inputs 1589b072555Sjeremylt const CeedScalar *X = in[0], *q_data = in[1]; 159ed264d09SValeria Barra // Outputs 160ed264d09SValeria Barra CeedScalar *true_soln = out[0], *rhs = out[1]; 161ed264d09SValeria Barra 162ed264d09SValeria Barra // Context 163ed264d09SValeria Barra const CeedScalar *context = (const CeedScalar *)ctx; 164ed264d09SValeria Barra const CeedScalar R = context[0]; 165ed264d09SValeria Barra 166ed264d09SValeria Barra // Quadrature Point Loop 167*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 168ed264d09SValeria Barra // Read global Cartesian coordinates 169ed264d09SValeria Barra CeedScalar x = X[i + Q * 0], y = X[i + Q * 1], z = X[i + Q * 2]; 170ed264d09SValeria Barra // Normalize quadrature point coordinates to sphere 171ed264d09SValeria Barra CeedScalar rad = sqrt(x * x + y * y + z * z); 172ed264d09SValeria Barra x *= R / rad; 173ed264d09SValeria Barra y *= R / rad; 174ed264d09SValeria Barra z *= R / rad; 175ed264d09SValeria Barra // Compute latitude and longitude 176ed264d09SValeria Barra const CeedScalar theta = asin(z / R); // latitude 177ed264d09SValeria Barra const CeedScalar lambda = atan2(y, x); // longitude 178ed264d09SValeria Barra 179ed264d09SValeria Barra true_soln[i + Q * 0] = sin(lambda) * cos(theta); 180ed264d09SValeria Barra 1819b072555Sjeremylt rhs[i + Q * 0] = q_data[i + Q * 0] * 2 * sin(lambda) * cos(theta) / (R * R); 182ed264d09SValeria Barra } // End of Quadrature Point Loop 183ed264d09SValeria Barra 184ed264d09SValeria Barra return 0; 185ed264d09SValeria Barra } 186ed264d09SValeria Barra 187e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 188ed264d09SValeria Barra // This QFunction applies the diffusion operator for a scalar field. 189ed264d09SValeria Barra // 190ed264d09SValeria Barra // Inputs: 191ed264d09SValeria Barra // ug - Input vector gradient at quadrature points 1929b072555Sjeremylt // q_data - Geometric factors 193ed264d09SValeria Barra // 194ed264d09SValeria Barra // Output: 195ed264d09SValeria Barra // vg - Output vector (test functions) gradient at quadrature points 196ed264d09SValeria Barra // 197ed264d09SValeria Barra // ----------------------------------------------------------------------------- 198*2b730f8bSJeremy L Thompson CEED_QFUNCTION(Diff)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 199ed264d09SValeria Barra // Inputs 2009b072555Sjeremylt const CeedScalar *ug = in[0], *q_data = in[1]; 201ed264d09SValeria Barra // Outputs 202ed264d09SValeria Barra CeedScalar *vg = out[0]; 203ed264d09SValeria Barra 204ed264d09SValeria Barra // Quadrature Point Loop 205*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 206ed264d09SValeria Barra // Read spatial derivatives of u 207*2b730f8bSJeremy L Thompson const CeedScalar du[2] = {ug[i + Q * 0], ug[i + Q * 1]}; 2089b072555Sjeremylt // Read q_data 2099b072555Sjeremylt const CeedScalar w_det_J = q_data[i + Q * 0]; 2109b072555Sjeremylt // -- Grad-to-Grad q_data 211ed264d09SValeria Barra // ---- dXdx_j,k * dXdx_k,j 212*2b730f8bSJeremy L Thompson const CeedScalar dXdxdXdx_T[2][2] = { 213*2b730f8bSJeremy L Thompson {q_data[i + Q * 1], q_data[i + Q * 3]}, 214*2b730f8bSJeremy L Thompson {q_data[i + Q * 3], q_data[i + Q * 2]} 215ed264d09SValeria Barra }; 216ed264d09SValeria Barra 217*2b730f8bSJeremy L Thompson for (int j = 0; j < 2; j++) { // j = direction of vg 218*2b730f8bSJeremy L Thompson vg[i + j * Q] = w_det_J * (du[0] * dXdxdXdx_T[0][j] + du[1] * dXdxdXdx_T[1][j]); 219*2b730f8bSJeremy L Thompson } 220ed264d09SValeria Barra } // End of Quadrature Point Loop 221ed264d09SValeria Barra 222ed264d09SValeria Barra return 0; 223ed264d09SValeria Barra } 224ed264d09SValeria Barra // ----------------------------------------------------------------------------- 225f6b55d2cSvaleriabarra 226f6b55d2cSvaleriabarra #endif // bp3sphere_h 227