xref: /honee/qfunctions/mass.h (revision 5a275bb67bb8312ddcb117417b4aaae9b609ecba)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Mass operator for Navier-Stokes example using PETSc
6 #include <ceed.h>
7 #include <math.h>
8 
9 // *****************************************************************************
10 // This QFunction applies the mass matrix to N components
11 //
12 // Inputs:
13 //   u      - Input vector at quadrature points
14 //   q_data - Quadrature weights
15 //
16 // Output:
17 //   v - Output vector at quadrature points
18 //
19 // *****************************************************************************
20 CEED_QFUNCTION_HELPER int Mass_N(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, const CeedInt N) {
21   const CeedScalar(*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
22   const CeedScalar(*q_data)        = in[1];
23   CeedScalar(*v)[CEED_Q_VLA]       = (CeedScalar(*)[CEED_Q_VLA])out[0];
24 
25   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
26     CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) { v[j][i] = q_data[i] * u[j][i]; }
27   }
28   return 0;
29 }
30 
31 CEED_QFUNCTION(Mass_1)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 1); }
32 
33 CEED_QFUNCTION(Mass_2)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 2); }
34 
35 CEED_QFUNCTION(Mass_3)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 3); }
36 
37 CEED_QFUNCTION(Mass_4)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 4); }
38 
39 CEED_QFUNCTION(Mass_5)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 5); }
40 
41 CEED_QFUNCTION(Mass_7)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 7); }
42 
43 CEED_QFUNCTION(Mass_9)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 9); }
44 
45 CEED_QFUNCTION(Mass_12)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 12); }
46 
47 CEED_QFUNCTION(Mass_22)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 22); }
48