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