1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2437c7c90SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3437c7c90SJeremy L Thompson //
4437c7c90SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5437c7c90SJeremy L Thompson //
6437c7c90SJeremy L Thompson // This file is part of CEED: http://github.com/ceed
7437c7c90SJeremy L Thompson
8c0b5abf0SJeremy L Thompson #include <ceed/types.h>
9437c7c90SJeremy L Thompson
apply(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)10437c7c90SJeremy L Thompson CEED_QFUNCTION(apply)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
11437c7c90SJeremy L Thompson // in[0] is gradient u_0, shape [2, num_comp=2, Q]
12437c7c90SJeremy L Thompson // in[1] is mass quadrature data, size (Q)
13437c7c90SJeremy L Thompson // in[2] is Poisson quadrature data, size (Q)
14437c7c90SJeremy L Thompson // in[3] is u_0, size (Q)
15437c7c90SJeremy L Thompson // in[4] is u_1, size (Q)
16437c7c90SJeremy L Thompson const CeedScalar *du_0 = in[0], *qd_mass = in[1], *qd_diff = in[2], *u_0 = in[3], *u_1 = in[4];
17437c7c90SJeremy L Thompson
18437c7c90SJeremy L Thompson // out[0] is output to multiply against v_0, size (Q)
19437c7c90SJeremy L Thompson // out[1] is output to multiply against v_1, size (Q)
20437c7c90SJeremy L Thompson // out[2] is output to multiply against gradient v_0, shape [2, num_comp=2, Q]
21437c7c90SJeremy L Thompson CeedScalar *v_0 = out[0], *v_1 = out[1], *dv_0 = out[2];
22437c7c90SJeremy L Thompson
23437c7c90SJeremy L Thompson const CeedInt num_comp_0 = 2;
24437c7c90SJeremy L Thompson
25437c7c90SJeremy L Thompson // Quadrature point loop
26437c7c90SJeremy L Thompson for (CeedInt i = 0; i < Q; i++) {
27437c7c90SJeremy L Thompson for (CeedInt c = 0; c < num_comp_0; c++) {
28437c7c90SJeremy L Thompson // Mass
29437c7c90SJeremy L Thompson v_0[i + Q * c] = qd_mass[i] * (c + 1) * u_0[i + Q * c];
30437c7c90SJeremy L Thompson // Diff
31437c7c90SJeremy L Thompson dv_0[i + Q * (0 * num_comp_0 + c)] =
32437c7c90SJeremy L Thompson qd_diff[i + Q * 0] * (c + 1) * du_0[i + Q * (0 * num_comp_0 + c)] + qd_diff[i + Q * 2] * du_0[i + Q * (1 * num_comp_0 + c)];
33437c7c90SJeremy L Thompson dv_0[i + Q * (1 * num_comp_0 + c)] =
34437c7c90SJeremy L Thompson qd_diff[i + Q * 2] * (c + 1) * du_0[i + Q * (0 * num_comp_0 + c)] + qd_diff[i + Q * 1] * du_0[i + Q * (1 * num_comp_0 + c)];
35437c7c90SJeremy L Thompson }
36437c7c90SJeremy L Thompson // Mass
37437c7c90SJeremy L Thompson v_1[i] = qd_mass[i] * u_1[i];
38437c7c90SJeremy L Thompson }
39437c7c90SJeremy L Thompson
40437c7c90SJeremy L Thompson return 0;
41437c7c90SJeremy L Thompson }
42