xref: /libCEED/examples/petsc/qfunctions/bps/bp3sphere.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 //
6ea61e9acSJeremy 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 
11c0b5abf0SJeremy L Thompson #include <ceed/types.h>
12c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS
13ed264d09SValeria Barra #include <math.h>
14c0b5abf0SJeremy L Thompson #endif
15ed264d09SValeria Barra 
16e83e87a5Sjeremylt // -----------------------------------------------------------------------------
17ea61e9acSJeremy L Thompson // This QFunction sets up the geometric factors required for integration and coordinate transformations when reference coordinates have a different
18ed264d09SValeria Barra // dimension than the one of physical coordinates
19ed264d09SValeria Barra //
20ed264d09SValeria Barra // Reference (parent) 2D coordinates: X \in [-1, 1]^2
21ed264d09SValeria Barra //
22ea61e9acSJeremy L Thompson // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 with R radius of the sphere
23ed264d09SValeria Barra //
24ea61e9acSJeremy 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
25ed264d09SValeria Barra //
26ed264d09SValeria Barra // Change of coordinates matrix computed by the library:
27ed264d09SValeria Barra //   (physical 3D coords relative to reference 2D coords)
28ed264d09SValeria Barra //   dxx_j/dX_i (indicial notation) [3 * 2]
29ed264d09SValeria Barra //
30ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D):
31ed264d09SValeria Barra //   dx_i/dxx_j (indicial notation) [3 * 3]
32ed264d09SValeria Barra //
33ed264d09SValeria Barra // Change of coordinates x (on the 2D manifold) relative to X (reference 2D):
34ed264d09SValeria Barra //   (by chain rule)
35ed264d09SValeria Barra //   dx_i/dX_j [3 * 2] = dx_i/dxx_k [3 * 3] * dxx_k/dX_j [3 * 2]
36ed264d09SValeria Barra //
379b072555Sjeremylt // mod_J is given by the magnitude of the cross product of the columns of dx_i/dX_j
38ed264d09SValeria Barra //
399b072555Sjeremylt // The quadrature data is stored in the array q_data.
40ed264d09SValeria Barra //
41ea61e9acSJeremy L Thompson // We require the determinant of the Jacobian to properly compute integrals of the form: int( u v )
42ed264d09SValeria Barra //
439b072555Sjeremylt // q_data[0]: mod_J * w
44ed264d09SValeria Barra //
45ea61e9acSJeremy 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
46ea61e9acSJeremy L Thompson // form: int( gradv gradu )
47ed264d09SValeria Barra //
48ed264d09SValeria Barra // dX_i/dx_j [2 * 3] = (dx_i/dX_j)+ = (dxdX^T dxdX)^(-1) dxdX
49ed264d09SValeria Barra //
50ac4340cfSJed Brown // and the product simplifies to yield the contravariant metric tensor
51ac4340cfSJed Brown //
52ac4340cfSJed Brown // g^{ij} = dX_i/dx_k dX_j/dx_k = (dxdX^T dxdX)^{-1}
53ac4340cfSJed Brown //
5408fade8cSvaleriabarra // Stored: g^{ij} (in Voigt convention) in
5508fade8cSvaleriabarra //
569b072555Sjeremylt //   q_data[1:3]: [dXdxdXdxT00 dXdxdXdxT01]
5708fade8cSvaleriabarra //               [dXdxdXdxT01 dXdxdXdxT11]
58ed264d09SValeria Barra // -----------------------------------------------------------------------------
SetupDiffGeo(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)592b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffGeo)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
60ed264d09SValeria Barra   const CeedScalar *X = in[0], *J = in[1], *w = in[2];
619b072555Sjeremylt   CeedScalar       *q_data = out[0];
62ed264d09SValeria Barra 
63ed264d09SValeria Barra   // Quadrature Point Loop
642b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
65ed264d09SValeria Barra     // Read global Cartesian coordinates
662b730f8bSJeremy L Thompson     const CeedScalar xx[3] = {X[i + 0 * Q], X[i + 1 * Q], X[i + 2 * Q]};
67ed264d09SValeria Barra 
68ed264d09SValeria Barra     // Read dxxdX Jacobian entries, stored as
69ed264d09SValeria Barra     // 0 3
70ed264d09SValeria Barra     // 1 4
71ed264d09SValeria Barra     // 2 5
722b730f8bSJeremy L Thompson     const CeedScalar dxxdX[3][2] = {
732b730f8bSJeremy L Thompson         {J[i + Q * 0], J[i + Q * 3]},
742b730f8bSJeremy L Thompson         {J[i + Q * 1], J[i + Q * 4]},
752b730f8bSJeremy L Thompson         {J[i + Q * 2], J[i + Q * 5]}
76ed264d09SValeria Barra     };
77ed264d09SValeria Barra 
78ed264d09SValeria Barra     // Setup
79ed264d09SValeria Barra     // x = xx (xx^T xx)^{-1/2}
80ed264d09SValeria Barra     // dx/dxx = I (xx^T xx)^{-1/2} - xx xx^T (xx^T xx)^{-3/2}
819b072555Sjeremylt     const CeedScalar mod_xx_sq = xx[0] * xx[0] + xx[1] * xx[1] + xx[2] * xx[2];
829b072555Sjeremylt     CeedScalar       xx_sq[3][3];
832b730f8bSJeremy L Thompson     for (int j = 0; j < 3; j++) {
842b730f8bSJeremy L Thompson       for (int k = 0; k < 3; k++) xx_sq[j][k] = xx[j] * xx[k] / (sqrt(mod_xx_sq) * mod_xx_sq);
852b730f8bSJeremy L Thompson     }
86ed264d09SValeria Barra 
872b730f8bSJeremy L Thompson     const CeedScalar dxdxx[3][3] = {
882b730f8bSJeremy L Thompson         {1. / sqrt(mod_xx_sq) - xx_sq[0][0], -xx_sq[0][1],                       -xx_sq[0][2]                      },
892b730f8bSJeremy L Thompson         {-xx_sq[1][0],                       1. / sqrt(mod_xx_sq) - xx_sq[1][1], -xx_sq[1][2]                      },
902b730f8bSJeremy L Thompson         {-xx_sq[2][0],                       -xx_sq[2][1],                       1. / sqrt(mod_xx_sq) - xx_sq[2][2]}
91ed264d09SValeria Barra     };
92ed264d09SValeria Barra 
93ed264d09SValeria Barra     CeedScalar dxdX[3][2];
942b730f8bSJeremy L Thompson     for (int j = 0; j < 3; j++) {
95ed264d09SValeria Barra       for (int k = 0; k < 2; k++) {
96ed264d09SValeria Barra         dxdX[j][k] = 0;
972b730f8bSJeremy L Thompson         for (int l = 0; l < 3; l++) dxdX[j][k] += dxdxx[j][l] * dxxdX[l][k];
982b730f8bSJeremy L Thompson       }
99ed264d09SValeria Barra     }
100ed264d09SValeria Barra 
101ed264d09SValeria Barra     // J is given by the cross product of the columns of dxdX
1022b730f8bSJeremy 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],
1032b730f8bSJeremy L Thompson                              dxdX[0][0] * dxdX[1][1] - dxdX[1][0] * dxdX[0][1]};
104ed264d09SValeria Barra 
105ed264d09SValeria Barra     // Use the magnitude of J as our detJ (volume scaling factor)
1069b072555Sjeremylt     const CeedScalar mod_J = sqrt(J[0] * J[0] + J[1] * J[1] + J[2] * J[2]);
107ed264d09SValeria Barra 
1089b072555Sjeremylt     // Interp-to-Interp q_data
1099b072555Sjeremylt     q_data[i + Q * 0] = mod_J * w[i];
110ed264d09SValeria Barra 
11108fade8cSvaleriabarra     // dxdX_k,j * dxdX_j,k
112ed264d09SValeria Barra     CeedScalar dxdXTdxdX[2][2];
1132b730f8bSJeremy L Thompson     for (int j = 0; j < 2; j++) {
114ed264d09SValeria Barra       for (int k = 0; k < 2; k++) {
115ed264d09SValeria Barra         dxdXTdxdX[j][k] = 0;
1162b730f8bSJeremy L Thompson         for (int l = 0; l < 3; l++) dxdXTdxdX[j][k] += dxdX[l][j] * dxdX[l][k];
1172b730f8bSJeremy L Thompson       }
118ed264d09SValeria Barra     }
119ed264d09SValeria Barra 
1202b730f8bSJeremy L Thompson     const CeedScalar detdxdXTdxdX = dxdXTdxdX[0][0] * dxdXTdxdX[1][1] - dxdXTdxdX[1][0] * dxdXTdxdX[0][1];
121ed264d09SValeria Barra 
12208fade8cSvaleriabarra     // Compute inverse of dxdXTdxdX, which is the 2x2 contravariant metric tensor g^{ij}
1239b072555Sjeremylt     CeedScalar dxdXTdxdX_inv[2][2];
1249b072555Sjeremylt     dxdXTdxdX_inv[0][0] = dxdXTdxdX[1][1] / detdxdXTdxdX;
1259b072555Sjeremylt     dxdXTdxdX_inv[0][1] = -dxdXTdxdX[0][1] / detdxdXTdxdX;
1269b072555Sjeremylt     dxdXTdxdX_inv[1][0] = -dxdXTdxdX[1][0] / detdxdXTdxdX;
1279b072555Sjeremylt     dxdXTdxdX_inv[1][1] = dxdXTdxdX[0][0] / detdxdXTdxdX;
128ed264d09SValeria Barra 
129ed264d09SValeria Barra     // Stored in Voigt convention
1309b072555Sjeremylt     q_data[i + Q * 1] = dxdXTdxdX_inv[0][0];
1319b072555Sjeremylt     q_data[i + Q * 2] = dxdXTdxdX_inv[1][1];
1329b072555Sjeremylt     q_data[i + Q * 3] = dxdXTdxdX_inv[0][1];
133ed264d09SValeria Barra   }  // End of Quadrature Point Loop
134ed264d09SValeria Barra 
135ed264d09SValeria Barra   // Return
136ed264d09SValeria Barra   return 0;
137ed264d09SValeria Barra }
138ed264d09SValeria Barra 
139e83e87a5Sjeremylt // -----------------------------------------------------------------------------
140ed264d09SValeria Barra // This QFunction sets up the rhs and true solution for the problem
141ed264d09SValeria Barra // -----------------------------------------------------------------------------
SetupDiffRhs(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)1422b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffRhs)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
143ed264d09SValeria Barra   // Inputs
1449b072555Sjeremylt   const CeedScalar *X = in[0], *q_data = in[1];
145ed264d09SValeria Barra   // Outputs
146ed264d09SValeria Barra   CeedScalar *true_soln = out[0], *rhs = out[1];
147ed264d09SValeria Barra 
148ed264d09SValeria Barra   // Context
149ed264d09SValeria Barra   const CeedScalar *context = (const CeedScalar *)ctx;
150ed264d09SValeria Barra   const CeedScalar  R       = context[0];
151ed264d09SValeria Barra 
152ed264d09SValeria Barra   // Quadrature Point Loop
1532b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
154ed264d09SValeria Barra     // Read global Cartesian coordinates
155ed264d09SValeria Barra     CeedScalar x = X[i + Q * 0], y = X[i + Q * 1], z = X[i + Q * 2];
156ed264d09SValeria Barra     // Normalize quadrature point coordinates to sphere
157ed264d09SValeria Barra     CeedScalar rad = sqrt(x * x + y * y + z * z);
158ed264d09SValeria Barra     x *= R / rad;
159ed264d09SValeria Barra     y *= R / rad;
160ed264d09SValeria Barra     z *= R / rad;
161ed264d09SValeria Barra     // Compute latitude and longitude
162ed264d09SValeria Barra     const CeedScalar theta  = asin(z / R);  // latitude
163ed264d09SValeria Barra     const CeedScalar lambda = atan2(y, x);  // longitude
164ed264d09SValeria Barra 
165ed264d09SValeria Barra     true_soln[i + Q * 0] = sin(lambda) * cos(theta);
166ed264d09SValeria Barra 
1679b072555Sjeremylt     rhs[i + Q * 0] = q_data[i + Q * 0] * 2 * sin(lambda) * cos(theta) / (R * R);
168ed264d09SValeria Barra   }  // End of Quadrature Point Loop
169ed264d09SValeria Barra 
170ed264d09SValeria Barra   return 0;
171ed264d09SValeria Barra }
172ed264d09SValeria Barra 
173e83e87a5Sjeremylt // -----------------------------------------------------------------------------
174ed264d09SValeria Barra // This QFunction applies the diffusion operator for a scalar field.
175ed264d09SValeria Barra //
176ed264d09SValeria Barra // Inputs:
177ed264d09SValeria Barra //   ug      - Input vector gradient at quadrature points
1789b072555Sjeremylt //   q_data  - Geometric factors
179ed264d09SValeria Barra //
180ed264d09SValeria Barra // Output:
181ed264d09SValeria Barra //   vg     - Output vector (test functions) gradient at quadrature points
182ed264d09SValeria Barra // -----------------------------------------------------------------------------
Diff(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)1832b730f8bSJeremy L Thompson CEED_QFUNCTION(Diff)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
184ed264d09SValeria Barra   // Inputs
1859b072555Sjeremylt   const CeedScalar *ug = in[0], *q_data = in[1];
186ed264d09SValeria Barra   // Outputs
187ed264d09SValeria Barra   CeedScalar *vg = out[0];
188ed264d09SValeria Barra 
189ed264d09SValeria Barra   // Quadrature Point Loop
1902b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
191ed264d09SValeria Barra     // Read spatial derivatives of u
1922b730f8bSJeremy L Thompson     const CeedScalar du[2] = {ug[i + Q * 0], ug[i + Q * 1]};
1939b072555Sjeremylt     // Read q_data
1949b072555Sjeremylt     const CeedScalar w_det_J = q_data[i + Q * 0];
1959b072555Sjeremylt     // -- Grad-to-Grad q_data
196ed264d09SValeria Barra     // ---- dXdx_j,k * dXdx_k,j
1972b730f8bSJeremy L Thompson     const CeedScalar dXdxdXdx_T[2][2] = {
1982b730f8bSJeremy L Thompson         {q_data[i + Q * 1], q_data[i + Q * 3]},
1992b730f8bSJeremy L Thompson         {q_data[i + Q * 3], q_data[i + Q * 2]}
200ed264d09SValeria Barra     };
201ed264d09SValeria Barra 
2022b730f8bSJeremy L Thompson     for (int j = 0; j < 2; j++) {  // j = direction of vg
2032b730f8bSJeremy L Thompson       vg[i + j * Q] = w_det_J * (du[0] * dXdxdXdx_T[0][j] + du[1] * dXdxdXdx_T[1][j]);
2042b730f8bSJeremy L Thompson     }
205ed264d09SValeria Barra   }  // End of Quadrature Point Loop
206ed264d09SValeria Barra 
207ed264d09SValeria Barra   return 0;
208ed264d09SValeria Barra }
209ed264d09SValeria Barra // -----------------------------------------------------------------------------
210