1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
10c0b5abf0SJeremy L Thompson #include <ceed/types.h>
11c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS
12c9c2c079SJeremy L Thompson #include <math.h>
13c0b5abf0SJeremy L Thompson #endif
1477841947SLeila Ghaffari
1577841947SLeila Ghaffari // *****************************************************************************
1677841947SLeila Ghaffari // This QFunction applies the mass matrix to five interlaced fields.
1777841947SLeila Ghaffari //
1877841947SLeila Ghaffari // Inputs:
1977841947SLeila Ghaffari // u - Input vector at quadrature points
2077841947SLeila Ghaffari // q_data - Quadrature weights
2177841947SLeila Ghaffari //
2277841947SLeila Ghaffari // Output:
2377841947SLeila Ghaffari // v - Output vector at quadrature points
2477841947SLeila Ghaffari //
2577841947SLeila Ghaffari // *****************************************************************************
Mass_N(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,const CeedInt N)26ef080ff9SJames Wright CEED_QFUNCTION_HELPER int Mass_N(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, const CeedInt N) {
2746603fc5SJames Wright const CeedScalar(*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2846603fc5SJames Wright const CeedScalar(*q_data) = in[1];
2977841947SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
3077841947SLeila Ghaffari
312b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
32ef080ff9SJames Wright CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) { v[j][i] = q_data[i] * u[j][i]; }
3377841947SLeila Ghaffari }
3477841947SLeila Ghaffari return 0;
3577841947SLeila Ghaffari }
3677841947SLeila Ghaffari
Mass_1(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)37ef080ff9SJames Wright CEED_QFUNCTION(Mass_1)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 1); }
38ef080ff9SJames Wright
Mass_5(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)39ef080ff9SJames Wright CEED_QFUNCTION(Mass_5)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 5); }
40ef080ff9SJames Wright
Mass_7(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)41052409adSJames Wright CEED_QFUNCTION(Mass_7)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 7); }
42052409adSJames Wright
Mass_9(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)43ef080ff9SJames Wright CEED_QFUNCTION(Mass_9)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 9); }
44ef080ff9SJames Wright
Mass_22(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)45ef080ff9SJames Wright CEED_QFUNCTION(Mass_22)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return Mass_N(ctx, Q, in, out, 22); }
46