xref: /libCEED/examples/petsc/qfunctions/area/areacube.h (revision 2459f3f1cd4d7d2e210e1c26d669bd2fde41a0b6)
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
7 
8 /// @file
9 /// libCEED QFunctions for mass operator example for a scalar field on the sphere using PETSc
10 
11 #ifndef areacube_h
12 #define areacube_h
13 
14 #include <math.h>
15 
16 // -----------------------------------------------------------------------------
17 // This QFunction sets up the geometric factor required for integration when
18 //   reference coordinates have a different dimension than the one of
19 //   physical coordinates
20 //
21 // Reference (parent) 2D coordinates: X \in [-1, 1]^2
22 //
23 // Global physical coordinates given by the mesh (3D): xx \in [-l, l]^3
24 //
25 // Local physical coordinates on the manifold (2D): x \in [-l, l]^2
26 //
27 // Change of coordinates matrix computed by the library:
28 //   (physical 3D coords relative to reference 2D coords)
29 //   dxx_j/dX_i (indicial notation) [3 * 2]
30 //
31 // Change of coordinates x (physical 2D) relative to xx (phyisical 3D):
32 //   dx_i/dxx_j (indicial notation) [2 * 3]
33 //
34 // Change of coordinates x (physical 2D) relative to X (reference 2D):
35 //   (by chain rule)
36 //   dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j
37 //
38 // The quadrature data is stored in the array q_data.
39 //
40 // We require the determinant of the Jacobian to properly compute integrals of
41 //   the form: int( u v )
42 //
43 // Qdata: w * det(dx_i/dX_j)
44 //
45 // -----------------------------------------------------------------------------
46 CEED_QFUNCTION(SetupMassGeoCube)(void *ctx, const CeedInt Q,
47                              const CeedScalar *const *in,
48                              CeedScalar *const *out) {
49   // Inputs
50   const CeedScalar *J = in[1], *w = in[2];
51   // Outputs
52   CeedScalar *q_data = out[0];
53 
54   // Quadrature Point Loop
55   CeedPragmaSIMD
56   for (CeedInt i=0; i<Q; i++) {
57     // Read dxxdX Jacobian entries, stored as
58     // 0 3
59     // 1 4
60     // 2 5
61     const CeedScalar dxxdX[3][2] = {{J[i+Q*0],
62                                      J[i+Q*3]},
63                                     {J[i+Q*1],
64                                      J[i+Q*4]},
65                                     {J[i+Q*2],
66                                      J[i+Q*5]}
67                                    };
68 
69     // Modulus of dxxdX column vectors
70     const CeedScalar mod_g_1 = sqrt(dxxdX[0][0]*dxxdX[0][0] +
71                                     dxxdX[1][0]*dxxdX[1][0] +
72                                     dxxdX[2][0]*dxxdX[2][0]);
73     const CeedScalar mod_g_2 = sqrt(dxxdX[0][1]*dxxdX[0][1] +
74                                     dxxdX[1][1]*dxxdX[1][1] +
75                                     dxxdX[2][1]*dxxdX[2][1]);
76 
77     // Use normalized column vectors of dxxdX as rows of dxdxx
78     const CeedScalar dxdxx[2][3] = {{dxxdX[0][0] / mod_g_1,
79                                      dxxdX[1][0] / mod_g_1,
80                                      dxxdX[2][0] / mod_g_1},
81                                     {dxxdX[0][1] / mod_g_2,
82                                      dxxdX[1][1] / mod_g_2,
83                                      dxxdX[2][1] / mod_g_2}
84                                    };
85 
86     CeedScalar dxdX[2][2];
87     for (int j=0; j<2; j++)
88       for (int k=0; k<2; k++) {
89         dxdX[j][k] = 0;
90         for (int l=0; l<3; l++)
91           dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k];
92       }
93 
94     q_data[i+Q*0] = (dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1]) * w[i];
95 
96   } // End of Quadrature Point Loop
97   return 0;
98 }
99 
100 // -----------------------------------------------------------------------------
101 // This QFunction applies the mass operator for a scalar field.
102 //
103 // Inputs:
104 //   u     - Input vector at quadrature points
105 //   q_data - Geometric factors
106 //
107 // Output:
108 //   v     - Output vector (test function) at quadrature points
109 //
110 // -----------------------------------------------------------------------------
111 CEED_QFUNCTION(Mass)(void *ctx, const CeedInt Q,
112                      const CeedScalar *const *in, CeedScalar *const *out) {
113   // Inputs
114   const CeedScalar *u = in[0], *q_data = in[1];
115   // Outputs
116   CeedScalar *v = out[0];
117 
118   // Quadrature Point Loop
119   CeedPragmaSIMD
120   for (CeedInt i=0; i<Q; i++)
121     v[i] = q_data[i] * u[i];
122 
123   return 0;
124 }
125 // -----------------------------------------------------------------------------
126 
127 #endif // areacube_h
128