xref: /libCEED/examples/fluids/qfunctions/mass.h (revision b0d170e7bc2e5c930ee481a47eb73044935a48a4)
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 /// Mass operator for Navier-Stokes example using PETSc
10 
11 #ifndef mass_h
12 #define mass_h
13 
14 #include <ceed.h>
15 #include <math.h>
16 
17 // *****************************************************************************
18 // This QFunction applies the mass matrix to five interlaced fields.
19 //
20 // Inputs:
21 //   u      - Input vector at quadrature points
22 //   q_data - Quadrature weights
23 //
24 // Output:
25 //   v - Output vector at quadrature points
26 //
27 // *****************************************************************************
28 CEED_QFUNCTION_HELPER int Mass_N(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, const CeedInt N) {
29   // Inputs
30   const CeedScalar(*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
31   const CeedScalar(*q_data)        = in[1];
32 
33   // Outputs
34   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
35 
36   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
37     CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) { v[j][i] = q_data[i] * u[j][i]; }
38   }
39   return 0;
40 }
41 
42 CEED_QFUNCTION(Mass_1)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 1); }
43 
44 CEED_QFUNCTION(Mass_5)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 5); }
45 
46 CEED_QFUNCTION(Mass_9)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 9); }
47 
48 CEED_QFUNCTION(Mass_22)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 22); }
49 // *****************************************************************************
50 
51 #endif  // mass_h
52