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 // 6*ea61e9acSJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7ed264d09SValeria Barra 8ed264d09SValeria Barra /// @file 9ed264d09SValeria Barra /// libCEED QFunctions for diffusion operator example for a scalar field on the sphere using PETSc 10ed264d09SValeria Barra 11f6b55d2cSvaleriabarra #ifndef bp3sphere_h 12f6b55d2cSvaleriabarra #define bp3sphere_h 13f6b55d2cSvaleriabarra 14c9c2c079SJeremy L Thompson #include <ceed.h> 15ed264d09SValeria Barra #include <math.h> 16ed264d09SValeria Barra 17e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 18*ea61e9acSJeremy L Thompson // This QFunction sets up the geometric factors required for integration and coordinate transformations when reference coordinates have a different 19ed264d09SValeria Barra // dimension than the one of physical coordinates 20ed264d09SValeria Barra // 21ed264d09SValeria Barra // Reference (parent) 2D coordinates: X \in [-1, 1]^2 22ed264d09SValeria Barra // 23*ea61e9acSJeremy L Thompson // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 with R radius of the sphere 24ed264d09SValeria Barra // 25*ea61e9acSJeremy L Thompson // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3 with l half edge of the cube inscribed in the sphere 26ed264d09SValeria Barra // 27ed264d09SValeria Barra // Change of coordinates matrix computed by the library: 28ed264d09SValeria Barra // (physical 3D coords relative to reference 2D coords) 29ed264d09SValeria Barra // dxx_j/dX_i (indicial notation) [3 * 2] 30ed264d09SValeria Barra // 31ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D): 32ed264d09SValeria Barra // dx_i/dxx_j (indicial notation) [3 * 3] 33ed264d09SValeria Barra // 34ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to X (reference 2D): 35ed264d09SValeria Barra // (by chain rule) 36ed264d09SValeria Barra // dx_i/dX_j [3 * 2] = dx_i/dxx_k [3 * 3] * dxx_k/dX_j [3 * 2] 37ed264d09SValeria Barra // 389b072555Sjeremylt // mod_J is given by the magnitude of the cross product of the columns of dx_i/dX_j 39ed264d09SValeria Barra // 409b072555Sjeremylt // The quadrature data is stored in the array q_data. 41ed264d09SValeria Barra // 42*ea61e9acSJeremy L Thompson // We require the determinant of the Jacobian to properly compute integrals of the form: int( u v ) 43ed264d09SValeria Barra // 449b072555Sjeremylt // q_data[0]: mod_J * w 45ed264d09SValeria Barra // 46*ea61e9acSJeremy L Thompson // We use the Moore–Penrose (left) pseudoinverse of dx_i/dX_j, to compute dX_i/dx_j (and its transpose), needed to properly compute integrals of the 47*ea61e9acSJeremy L Thompson // form: int( gradv gradu ) 48ed264d09SValeria Barra // 49ed264d09SValeria Barra // dX_i/dx_j [2 * 3] = (dx_i/dX_j)+ = (dxdX^T dxdX)^(-1) dxdX 50ed264d09SValeria Barra // 51ac4340cfSJed Brown // and the product simplifies to yield the contravariant metric tensor 52ac4340cfSJed Brown // 53ac4340cfSJed Brown // g^{ij} = dX_i/dx_k dX_j/dx_k = (dxdX^T dxdX)^{-1} 54ac4340cfSJed Brown // 5508fade8cSvaleriabarra // Stored: g^{ij} (in Voigt convention) in 5608fade8cSvaleriabarra // 579b072555Sjeremylt // q_data[1:3]: [dXdxdXdxT00 dXdxdXdxT01] 5808fade8cSvaleriabarra // [dXdxdXdxT01 dXdxdXdxT11] 59ed264d09SValeria Barra // ----------------------------------------------------------------------------- 602b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffGeo)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 61ed264d09SValeria Barra const CeedScalar *X = in[0], *J = in[1], *w = in[2]; 629b072555Sjeremylt CeedScalar *q_data = out[0]; 63ed264d09SValeria Barra 64ed264d09SValeria Barra // Quadrature Point Loop 652b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 66ed264d09SValeria Barra // Read global Cartesian coordinates 672b730f8bSJeremy L Thompson const CeedScalar xx[3] = {X[i + 0 * Q], X[i + 1 * Q], X[i + 2 * Q]}; 68ed264d09SValeria Barra 69ed264d09SValeria Barra // Read dxxdX Jacobian entries, stored as 70ed264d09SValeria Barra // 0 3 71ed264d09SValeria Barra // 1 4 72ed264d09SValeria Barra // 2 5 732b730f8bSJeremy L Thompson const CeedScalar dxxdX[3][2] = { 742b730f8bSJeremy L Thompson {J[i + Q * 0], J[i + Q * 3]}, 752b730f8bSJeremy L Thompson {J[i + Q * 1], J[i + Q * 4]}, 762b730f8bSJeremy L Thompson {J[i + Q * 2], J[i + Q * 5]} 77ed264d09SValeria Barra }; 78ed264d09SValeria Barra 79ed264d09SValeria Barra // Setup 80ed264d09SValeria Barra // x = xx (xx^T xx)^{-1/2} 81ed264d09SValeria Barra // dx/dxx = I (xx^T xx)^{-1/2} - xx xx^T (xx^T xx)^{-3/2} 829b072555Sjeremylt const CeedScalar mod_xx_sq = xx[0] * xx[0] + xx[1] * xx[1] + xx[2] * xx[2]; 839b072555Sjeremylt CeedScalar xx_sq[3][3]; 842b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) { 852b730f8bSJeremy L Thompson for (int k = 0; k < 3; k++) xx_sq[j][k] = xx[j] * xx[k] / (sqrt(mod_xx_sq) * mod_xx_sq); 862b730f8bSJeremy L Thompson } 87ed264d09SValeria Barra 882b730f8bSJeremy L Thompson const CeedScalar dxdxx[3][3] = { 892b730f8bSJeremy L Thompson {1. / sqrt(mod_xx_sq) - xx_sq[0][0], -xx_sq[0][1], -xx_sq[0][2] }, 902b730f8bSJeremy L Thompson {-xx_sq[1][0], 1. / sqrt(mod_xx_sq) - xx_sq[1][1], -xx_sq[1][2] }, 912b730f8bSJeremy L Thompson {-xx_sq[2][0], -xx_sq[2][1], 1. / sqrt(mod_xx_sq) - xx_sq[2][2]} 92ed264d09SValeria Barra }; 93ed264d09SValeria Barra 94ed264d09SValeria Barra CeedScalar dxdX[3][2]; 952b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) { 96ed264d09SValeria Barra for (int k = 0; k < 2; k++) { 97ed264d09SValeria Barra dxdX[j][k] = 0; 982b730f8bSJeremy L Thompson for (int l = 0; l < 3; l++) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k]; 992b730f8bSJeremy L Thompson } 100ed264d09SValeria Barra } 101ed264d09SValeria Barra 102ed264d09SValeria Barra // J is given by the cross product of the columns of dxdX 1032b730f8bSJeremy 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], 1042b730f8bSJeremy L Thompson dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]}; 105ed264d09SValeria Barra 106ed264d09SValeria Barra // Use the magnitude of J as our detJ (volume scaling factor) 1079b072555Sjeremylt const CeedScalar mod_J = sqrt(J[0] * J[0] + J[1] * J[1] + J[2] * J[2]); 108ed264d09SValeria Barra 1099b072555Sjeremylt // Interp-to-Interp q_data 1109b072555Sjeremylt q_data[i + Q * 0] = mod_J * w[i]; 111ed264d09SValeria Barra 11208fade8cSvaleriabarra // dxdX_k,j * dxdX_j,k 113ed264d09SValeria Barra CeedScalar dxdXTdxdX[2][2]; 1142b730f8bSJeremy L Thompson for (int j = 0; j < 2; j++) { 115ed264d09SValeria Barra for (int k = 0; k < 2; k++) { 116ed264d09SValeria Barra dxdXTdxdX[j][k] = 0; 1172b730f8bSJeremy L Thompson for (int l = 0; l < 3; l++) dxdXTdxdX[j][k] += dxdX[l][j] * dxdX[l][k]; 1182b730f8bSJeremy L Thompson } 119ed264d09SValeria Barra } 120ed264d09SValeria Barra 1212b730f8bSJeremy L Thompson const CeedScalar detdxdXTdxdX = dxdXTdxdX[0][0] * dxdXTdxdX[1][1] - dxdXTdxdX[1][0] * dxdXTdxdX[0][1]; 122ed264d09SValeria Barra 12308fade8cSvaleriabarra // Compute inverse of dxdXTdxdX, which is the 2x2 contravariant metric tensor g^{ij} 1249b072555Sjeremylt CeedScalar dxdXTdxdX_inv[2][2]; 1259b072555Sjeremylt dxdXTdxdX_inv[0][0] = dxdXTdxdX[1][1] / detdxdXTdxdX; 1269b072555Sjeremylt dxdXTdxdX_inv[0][1] = -dxdXTdxdX[0][1] / detdxdXTdxdX; 1279b072555Sjeremylt dxdXTdxdX_inv[1][0] = -dxdXTdxdX[1][0] / detdxdXTdxdX; 1289b072555Sjeremylt dxdXTdxdX_inv[1][1] = dxdXTdxdX[0][0] / detdxdXTdxdX; 129ed264d09SValeria Barra 130ed264d09SValeria Barra // Stored in Voigt convention 1319b072555Sjeremylt q_data[i + Q * 1] = dxdXTdxdX_inv[0][0]; 1329b072555Sjeremylt q_data[i + Q * 2] = dxdXTdxdX_inv[1][1]; 1339b072555Sjeremylt q_data[i + Q * 3] = dxdXTdxdX_inv[0][1]; 134ed264d09SValeria Barra } // End of Quadrature Point Loop 135ed264d09SValeria Barra 136ed264d09SValeria Barra // Return 137ed264d09SValeria Barra return 0; 138ed264d09SValeria Barra } 139ed264d09SValeria Barra 140e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 141ed264d09SValeria Barra // This QFunction sets up the rhs and true solution for the problem 142ed264d09SValeria Barra // ----------------------------------------------------------------------------- 1432b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffRhs)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 144ed264d09SValeria Barra // Inputs 1459b072555Sjeremylt const CeedScalar *X = in[0], *q_data = in[1]; 146ed264d09SValeria Barra // Outputs 147ed264d09SValeria Barra CeedScalar *true_soln = out[0], *rhs = out[1]; 148ed264d09SValeria Barra 149ed264d09SValeria Barra // Context 150ed264d09SValeria Barra const CeedScalar *context = (const CeedScalar *)ctx; 151ed264d09SValeria Barra const CeedScalar R = context[0]; 152ed264d09SValeria Barra 153ed264d09SValeria Barra // Quadrature Point Loop 1542b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 155ed264d09SValeria Barra // Read global Cartesian coordinates 156ed264d09SValeria Barra CeedScalar x = X[i + Q * 0], y = X[i + Q * 1], z = X[i + Q * 2]; 157ed264d09SValeria Barra // Normalize quadrature point coordinates to sphere 158ed264d09SValeria Barra CeedScalar rad = sqrt(x * x + y * y + z * z); 159ed264d09SValeria Barra x *= R / rad; 160ed264d09SValeria Barra y *= R / rad; 161ed264d09SValeria Barra z *= R / rad; 162ed264d09SValeria Barra // Compute latitude and longitude 163ed264d09SValeria Barra const CeedScalar theta = asin(z / R); // latitude 164ed264d09SValeria Barra const CeedScalar lambda = atan2(y, x); // longitude 165ed264d09SValeria Barra 166ed264d09SValeria Barra true_soln[i + Q * 0] = sin(lambda) * cos(theta); 167ed264d09SValeria Barra 1689b072555Sjeremylt rhs[i + Q * 0] = q_data[i + Q * 0] * 2 * sin(lambda) * cos(theta) / (R * R); 169ed264d09SValeria Barra } // End of Quadrature Point Loop 170ed264d09SValeria Barra 171ed264d09SValeria Barra return 0; 172ed264d09SValeria Barra } 173ed264d09SValeria Barra 174e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 175ed264d09SValeria Barra // This QFunction applies the diffusion operator for a scalar field. 176ed264d09SValeria Barra // 177ed264d09SValeria Barra // Inputs: 178ed264d09SValeria Barra // ug - Input vector gradient at quadrature points 1799b072555Sjeremylt // q_data - Geometric factors 180ed264d09SValeria Barra // 181ed264d09SValeria Barra // Output: 182ed264d09SValeria Barra // vg - Output vector (test functions) gradient at quadrature points 183ed264d09SValeria Barra // ----------------------------------------------------------------------------- 1842b730f8bSJeremy L Thompson CEED_QFUNCTION(Diff)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 185ed264d09SValeria Barra // Inputs 1869b072555Sjeremylt const CeedScalar *ug = in[0], *q_data = in[1]; 187ed264d09SValeria Barra // Outputs 188ed264d09SValeria Barra CeedScalar *vg = out[0]; 189ed264d09SValeria Barra 190ed264d09SValeria Barra // Quadrature Point Loop 1912b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 192ed264d09SValeria Barra // Read spatial derivatives of u 1932b730f8bSJeremy L Thompson const CeedScalar du[2] = {ug[i + Q * 0], ug[i + Q * 1]}; 1949b072555Sjeremylt // Read q_data 1959b072555Sjeremylt const CeedScalar w_det_J = q_data[i + Q * 0]; 1969b072555Sjeremylt // -- Grad-to-Grad q_data 197ed264d09SValeria Barra // ---- dXdx_j,k * dXdx_k,j 1982b730f8bSJeremy L Thompson const CeedScalar dXdxdXdx_T[2][2] = { 1992b730f8bSJeremy L Thompson {q_data[i + Q * 1], q_data[i + Q * 3]}, 2002b730f8bSJeremy L Thompson {q_data[i + Q * 3], q_data[i + Q * 2]} 201ed264d09SValeria Barra }; 202ed264d09SValeria Barra 2032b730f8bSJeremy L Thompson for (int j = 0; j < 2; j++) { // j = direction of vg 2042b730f8bSJeremy L Thompson vg[i + j * Q] = w_det_J * (du[0] * dXdxdXdx_T[0][j] + du[1] * dXdxdXdx_T[1][j]); 2052b730f8bSJeremy L Thompson } 206ed264d09SValeria Barra } // End of Quadrature Point Loop 207ed264d09SValeria Barra 208ed264d09SValeria Barra return 0; 209ed264d09SValeria Barra } 210ed264d09SValeria Barra // ----------------------------------------------------------------------------- 211f6b55d2cSvaleriabarra 212f6b55d2cSvaleriabarra #endif // bp3sphere_h 213