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