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.
3cb32e2e7SValeria Barra //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5cb32e2e7SValeria Barra //
63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed
7cb32e2e7SValeria Barra
8cb32e2e7SValeria Barra /// @file
9cb32e2e7SValeria Barra /// libCEED QFunctions for mass operator example using PETSc
10cb32e2e7SValeria Barra
11c0b5abf0SJeremy L Thompson #include <ceed/types.h>
12c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS
13cb32e2e7SValeria Barra #include <math.h>
14c0b5abf0SJeremy L Thompson #endif
15cb32e2e7SValeria Barra
16e83e87a5Sjeremylt // -----------------------------------------------------------------------------
17ed264d09SValeria Barra // This QFunction sets up the rhs and true solution for the problem
18cb32e2e7SValeria Barra // -----------------------------------------------------------------------------
SetupMassRhs3(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)192b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupMassRhs3)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
20e83e87a5Sjeremylt const CeedScalar *x = in[0], *w = in[1];
21cb32e2e7SValeria Barra CeedScalar *true_soln = out[0], *rhs = out[1];
22cb32e2e7SValeria Barra
23cb32e2e7SValeria Barra // Quadrature Point Loop
242b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
25cb32e2e7SValeria Barra // Component 1
26cb32e2e7SValeria Barra true_soln[i + 0 * Q] = sqrt(x[i] * x[i] + x[i + Q] * x[i + Q] + x[i + 2 * Q] * x[i + 2 * Q]);
27cb32e2e7SValeria Barra // Component 2
2882311801Sjeremylt true_soln[i + 1 * Q] = 2 * true_soln[i + 0 * Q];
29cb32e2e7SValeria Barra // Component 3
3082311801Sjeremylt true_soln[i + 2 * Q] = 3 * true_soln[i + 0 * Q];
31cb32e2e7SValeria Barra
32cb32e2e7SValeria Barra // Component 1
33e83e87a5Sjeremylt rhs[i + 0 * Q] = w[i] * true_soln[i + 0 * Q];
34cb32e2e7SValeria Barra // Component 2
3582311801Sjeremylt rhs[i + 1 * Q] = 2 * rhs[i + 0 * Q];
36cb32e2e7SValeria Barra // Component 3
3782311801Sjeremylt rhs[i + 2 * Q] = 3 * rhs[i + 0 * Q];
38cb32e2e7SValeria Barra } // End of Quadrature Point Loop
39cb32e2e7SValeria Barra return 0;
40cb32e2e7SValeria Barra }
41cb32e2e7SValeria Barra
42e83e87a5Sjeremylt // -----------------------------------------------------------------------------
43ed264d09SValeria Barra // This QFunction applies the mass operator for a vector field of 3 components.
44ed264d09SValeria Barra //
45ed264d09SValeria Barra // Inputs:
46ed264d09SValeria Barra // u - Input vector at quadrature points
479b072555Sjeremylt // q_data - Geometric factors
48ed264d09SValeria Barra //
49ed264d09SValeria Barra // Output:
50ed264d09SValeria Barra // v - Output vector (test functions) at quadrature points
51cb32e2e7SValeria Barra // -----------------------------------------------------------------------------
Mass3(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)522b730f8bSJeremy L Thompson CEED_QFUNCTION(Mass3)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
539b072555Sjeremylt const CeedScalar *u = in[0], *q_data = in[1];
54cb32e2e7SValeria Barra CeedScalar *v = out[0];
55cb32e2e7SValeria Barra
56cb32e2e7SValeria Barra // Quadrature Point Loop
572b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
58cb32e2e7SValeria Barra // Component 1
599b072555Sjeremylt v[i + 0 * Q] = q_data[i] * u[i + 0 * Q];
60cb32e2e7SValeria Barra // Component 2
619b072555Sjeremylt v[i + 1 * Q] = q_data[i] * u[i + 1 * Q];
62cb32e2e7SValeria Barra // Component 3
639b072555Sjeremylt v[i + 2 * Q] = q_data[i] * u[i + 2 * Q];
64cb32e2e7SValeria Barra } // End of Quadrature Point Loop
65cb32e2e7SValeria Barra return 0;
66cb32e2e7SValeria Barra }
67cb32e2e7SValeria Barra // -----------------------------------------------------------------------------
68