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