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