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