xref: /libCEED/examples/fluids/qfunctions/mass.h (revision ef080ff9ce83a1650979d1b767b88f0d6e3ee6ca)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Mass operator for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #ifndef mass_h
1277841947SLeila Ghaffari #define mass_h
1377841947SLeila Ghaffari 
14ba6664aeSJames Wright #include <ceed.h>
15c9c2c079SJeremy L Thompson #include <math.h>
1677841947SLeila Ghaffari 
1777841947SLeila Ghaffari // *****************************************************************************
1877841947SLeila Ghaffari // This QFunction applies the mass matrix to five interlaced fields.
1977841947SLeila Ghaffari //
2077841947SLeila Ghaffari // Inputs:
2177841947SLeila Ghaffari //   u      - Input vector at quadrature points
2277841947SLeila Ghaffari //   q_data - Quadrature weights
2377841947SLeila Ghaffari //
2477841947SLeila Ghaffari // Output:
2577841947SLeila Ghaffari //   v - Output vector at quadrature points
2677841947SLeila Ghaffari //
2777841947SLeila Ghaffari // *****************************************************************************
28*ef080ff9SJames Wright CEED_QFUNCTION_HELPER int Mass_N(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, const CeedInt N) {
2977841947SLeila Ghaffari   // Inputs
3046603fc5SJames Wright   const CeedScalar(*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
3146603fc5SJames Wright   const CeedScalar(*q_data)        = in[1];
3277841947SLeila Ghaffari 
3377841947SLeila Ghaffari   // Outputs
3477841947SLeila Ghaffari   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
3577841947SLeila Ghaffari 
362b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
37*ef080ff9SJames Wright     CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) { v[j][i] = q_data[i] * u[j][i]; }
3877841947SLeila Ghaffari   }
3977841947SLeila Ghaffari   return 0;
4077841947SLeila Ghaffari }
4177841947SLeila Ghaffari 
42*ef080ff9SJames Wright CEED_QFUNCTION(Mass_1)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 1); }
43*ef080ff9SJames Wright 
44*ef080ff9SJames Wright CEED_QFUNCTION(Mass_5)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 5); }
45*ef080ff9SJames Wright 
46*ef080ff9SJames Wright CEED_QFUNCTION(Mass_9)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 9); }
47*ef080ff9SJames Wright 
48*ef080ff9SJames Wright CEED_QFUNCTION(Mass_22)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 22); }
4977841947SLeila Ghaffari // *****************************************************************************
5077841947SLeila Ghaffari 
5177841947SLeila Ghaffari #endif  // mass_h
52