xref: /libCEED/examples/python/qfunctions/ex1-volume.h (revision d89549c6e0afe634fc6cadb472a2a3224f4b1b63)
1 // Copyright (c) 2017-2025, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 #pragma once
8 
9 #include <ceed/types.h>
10 #include "ex-common.h"
11 
12 /// libCEED Q-function for building quadrature data for a mass operator
13 CEED_QFUNCTION(build_mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
14   // in[0] is Jacobians with shape [dim, dim, Q]
15   // in[1] is quadrature weights with shape [1, Q]
16   const CeedScalar    *w          = in[1];
17   CeedScalar          *q_data     = out[0];
18   struct BuildContext *build_data = (struct BuildContext *)ctx;
19 
20   switch (build_data->dim + 10 * build_data->space_dim) {
21     case 11: {
22       const CeedScalar(*J)[1][CEED_Q_VLA] = (const CeedScalar(*)[1][CEED_Q_VLA])in[0];
23 
24       // Quadrature Point Loop
25       CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { q_data[i] = J[0][0][i] * w[i]; }  // End of Quadrature Point Loop
26     } break;
27     case 22: {
28       const CeedScalar(*J)[2][CEED_Q_VLA] = (const CeedScalar(*)[2][CEED_Q_VLA])in[0];
29 
30       // Quadrature Point Loop
31       CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
32         q_data[i] = (J[0][0][i] * J[1][1][i] - J[0][1][i] * J[1][0][i]) * w[i];
33       }  // End of Quadrature Point Loop
34     } break;
35     case 33: {
36       const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0];
37 
38       // Quadrature Point Loop
39       CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
40         q_data[i] =
41             (J[0][0][i] * (J[1][1][i] * J[2][2][i] - J[1][2][i] * J[2][1][i]) - J[0][1][i] * (J[1][0][i] * J[2][2][i] - J[1][2][i] * J[2][0][i]) +
42              J[0][2][i] * (J[1][0][i] * J[2][1][i] - J[1][1][i] * J[2][0][i])) *
43             w[i];
44       }  // End of Quadrature Point Loop
45     } break;
46   }
47   return CEED_ERROR_SUCCESS;
48 }
49 
50 /// libCEED Q-function for applying a mass operator
51 CEED_QFUNCTION(apply_mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
52   // in[0], out[0] are solution variables with shape [1, Q]
53   // in[1] is quadrature data with shape [1, Q]
54   const CeedScalar *u = in[0], *q_data = in[1];
55   CeedScalar       *v = out[0];
56 
57   // Quadrature Point Loop
58   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { v[i] = q_data[i] * u[i]; }  // End of Quadrature Point Loop
59   return CEED_ERROR_SUCCESS;
60 }
61