Line data Source code
1 : // Copyright Amazon.com, Inc. or its affiliates. All Rights Reserved.
2 : // SPDX-License-Identifier: Apache-2.0
3 :
4 : #ifndef PALACE_LIBCEED_HCURL_MASS_33_QF_H
5 : #define PALACE_LIBCEED_HCURL_MASS_33_QF_H
6 :
7 : #include "../coeff/coeff_1_qf.h"
8 : #include "../coeff/coeff_3_qf.h"
9 : #include "utils_33_qf.h"
10 :
11 88 : CEED_QFUNCTION(f_apply_hcurlmass_33)(void *__restrict__ ctx, CeedInt Q,
12 : const CeedScalar *const *in, CeedScalar *const *out)
13 : {
14 88 : const CeedScalar *attr = in[0], *wdetJ = in[0] + Q, *adjJt = in[0] + 2 * Q, *u = in[1],
15 88 : *gradu = in[2];
16 88 : CeedScalar *__restrict__ v = out[0], *__restrict__ gradv = out[1];
17 :
18 46104 : CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++)
19 : {
20 : {
21 92032 : const CeedScalar coeff = CoeffUnpack1((const CeedIntScalar *)ctx, (CeedInt)attr[i]);
22 :
23 46016 : v[i] = coeff * wdetJ[i] * u[i];
24 : }
25 : {
26 46016 : const CeedScalar u_loc[3] = {gradu[i + Q * 0], gradu[i + Q * 1], gradu[i + Q * 2]};
27 : CeedScalar coeff[9], adjJt_loc[9], v_loc[3];
28 46016 : CoeffUnpack3(CoeffPairSecond<1>((const CeedIntScalar *)ctx), (CeedInt)attr[i], coeff);
29 46016 : MatUnpack33(adjJt + i, Q, adjJt_loc);
30 : MultAtBCx33(adjJt_loc, coeff, adjJt_loc, u_loc, v_loc);
31 :
32 46016 : gradv[i + Q * 0] = wdetJ[i] * v_loc[0];
33 46016 : gradv[i + Q * 1] = wdetJ[i] * v_loc[1];
34 46016 : gradv[i + Q * 2] = wdetJ[i] * v_loc[2];
35 : }
36 : }
37 88 : return 0;
38 : }
39 :
40 : #endif // PALACE_LIBCEED_HCURL_MASS_33_QF_H
|