1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2b6972d74SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3b6972d74SZach Atkins //
4b6972d74SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
5b6972d74SZach Atkins //
6b6972d74SZach Atkins // This file is part of CEED: http://github.com/ceed
798285ab4SZach Atkins
8c0b5abf0SJeremy L Thompson #include <ceed/types.h>
9b6972d74SZach Atkins
SetupMass(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)10b6972d74SZach Atkins CEED_QFUNCTION(SetupMass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
11b6972d74SZach Atkins const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0];
12b6972d74SZach Atkins const CeedScalar *w = in[1];
13b6972d74SZach Atkins CeedScalar *q_data = out[0];
14b6972d74SZach Atkins
15b6972d74SZach Atkins // Quadrature point loop
16b6972d74SZach Atkins CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
17b6972d74SZach Atkins const CeedScalar B11 = J[1][1][i] * J[2][2][i] - J[1][2][i] * J[2][1][i];
18b6972d74SZach Atkins const CeedScalar B12 = J[0][2][i] * J[2][1][i] - J[0][1][i] * J[2][2][i];
19b6972d74SZach Atkins const CeedScalar B13 = J[0][1][i] * J[1][2][i] - J[0][2][i] * J[1][1][i];
20b6972d74SZach Atkins
21b6972d74SZach Atkins q_data[i] = w[i] * (J[0][0][i] * B11 + J[1][0][i] * B12 + J[2][0][i] * B13);
22b6972d74SZach Atkins }
23b6972d74SZach Atkins return 0;
24b6972d74SZach Atkins }
25b6972d74SZach Atkins
Mass(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)26b6972d74SZach Atkins CEED_QFUNCTION(Mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
27b6972d74SZach Atkins const CeedScalar *rho = (const CeedScalar *)in[0];
28b6972d74SZach Atkins const CeedScalar(*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
29b6972d74SZach Atkins CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
30b6972d74SZach Atkins
31b6972d74SZach Atkins const CeedInt num_comp = *(CeedInt *)ctx;
32b6972d74SZach Atkins // Quadrature point loop
33b6972d74SZach Atkins CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
34b6972d74SZach Atkins for (CeedInt j = 0; j < num_comp; j++) v[j][i] = rho[i] * u[j][i];
35b6972d74SZach Atkins }
36b6972d74SZach Atkins return 0;
37b6972d74SZach Atkins }
38