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