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.
3ed264d09SValeria Barra //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5ed264d09SValeria Barra //
63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed
7ed264d09SValeria Barra
8ed264d09SValeria Barra /// @file
9ed264d09SValeria Barra /// libCEED QFunctions for mass operator example for a vector field on the sphere using PETSc
10ed264d09SValeria Barra
11c0b5abf0SJeremy L Thompson #include <ceed/types.h>
12c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS
13ed264d09SValeria Barra #include <math.h>
14c0b5abf0SJeremy L Thompson #endif
15ed264d09SValeria Barra
16e83e87a5Sjeremylt // -----------------------------------------------------------------------------
17ed264d09SValeria Barra // This QFunction sets up the rhs and true solution for the problem
18ed264d09SValeria Barra // -----------------------------------------------------------------------------
SetupDiffRhs3(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)192b730f8bSJeremy L Thompson CEED_QFUNCTION(SetupDiffRhs3)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
20ed264d09SValeria Barra // Inputs
219b072555Sjeremylt const CeedScalar *X = in[0], *q_data = in[1];
22ed264d09SValeria Barra // Outputs
23ed264d09SValeria Barra CeedScalar *true_soln = out[0], *rhs = out[1];
24ed264d09SValeria Barra
25ed264d09SValeria Barra // Context
26ed264d09SValeria Barra const CeedScalar *context = (const CeedScalar *)ctx;
27ed264d09SValeria Barra const CeedScalar R = context[0];
28ed264d09SValeria Barra
29ed264d09SValeria Barra // Quadrature Point Loop
302b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
31ed264d09SValeria Barra // Read global Cartesian coordinates
32ed264d09SValeria Barra CeedScalar x = X[i + Q * 0], y = X[i + Q * 1], z = X[i + Q * 2];
33ed264d09SValeria Barra // Normalize quadrature point coordinates to sphere
34ed264d09SValeria Barra CeedScalar rad = sqrt(x * x + y * y + z * z);
35ed264d09SValeria Barra x *= R / rad;
36ed264d09SValeria Barra y *= R / rad;
37ed264d09SValeria Barra z *= R / rad;
38ed264d09SValeria Barra // Compute latitude and longitude
39ed264d09SValeria Barra const CeedScalar theta = asin(z / R); // latitude
40ed264d09SValeria Barra const CeedScalar lambda = atan2(y, x); // longitude
41ed264d09SValeria Barra
429b072555Sjeremylt // Use absolute value of latitude for true solution
43ed264d09SValeria Barra // Component 1
44ed264d09SValeria Barra true_soln[i + 0 * Q] = sin(lambda) * cos(theta);
45ed264d09SValeria Barra // Component 2
46ed264d09SValeria Barra true_soln[i + 1 * Q] = 2 * true_soln[i + 0 * Q];
47ed264d09SValeria Barra // Component 3
48ed264d09SValeria Barra true_soln[i + 2 * Q] = 3 * true_soln[i + 0 * Q];
49ed264d09SValeria Barra
50ed264d09SValeria Barra // Component 1
519b072555Sjeremylt rhs[i + 0 * Q] = q_data[i + Q * 0] * 2 * sin(lambda) * cos(theta) / (R * R);
52ed264d09SValeria Barra // Component 2
53ed264d09SValeria Barra rhs[i + 1 * Q] = 2 * rhs[i + 0 * Q];
54ed264d09SValeria Barra // Component 3
55ed264d09SValeria Barra rhs[i + 2 * Q] = 3 * rhs[i + 0 * Q];
56ed264d09SValeria Barra } // End of Quadrature Point Loop
57ed264d09SValeria Barra
58ed264d09SValeria Barra return 0;
59ed264d09SValeria Barra }
60ed264d09SValeria Barra
61e83e87a5Sjeremylt // -----------------------------------------------------------------------------
62ed264d09SValeria Barra // This QFunction applies the diffusion operator for a vector field of 3 components.
63ed264d09SValeria Barra //
64ed264d09SValeria Barra // Inputs:
65ed264d09SValeria Barra // ug - Input vector Jacobian at quadrature points
669b072555Sjeremylt // q_data - Geometric factors
67ed264d09SValeria Barra //
68ed264d09SValeria Barra // Output:
69ed264d09SValeria Barra // vJ - Output vector (test functions) Jacobian at quadrature points
70ed264d09SValeria Barra // -----------------------------------------------------------------------------
Diff3(void * ctx,const CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)712b730f8bSJeremy L Thompson CEED_QFUNCTION(Diff3)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
729b072555Sjeremylt const CeedScalar *ug = in[0], *q_data = in[1];
73ed264d09SValeria Barra CeedScalar *vJ = out[0];
74ed264d09SValeria Barra
75ed264d09SValeria Barra // Quadrature Point Loop
762b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
77ed264d09SValeria Barra // Read spatial derivatives of u
782b730f8bSJeremy L Thompson const CeedScalar uJ[3][2] = {
792b730f8bSJeremy L Thompson {ug[i + (0 + 0 * 3) * Q], ug[i + (0 + 1 * 3) * Q]},
802b730f8bSJeremy L Thompson {ug[i + (1 + 0 * 3) * Q], ug[i + (1 + 1 * 3) * Q]},
812b730f8bSJeremy L Thompson {ug[i + (2 + 0 * 3) * Q], ug[i + (2 + 1 * 3) * Q]}
82ed264d09SValeria Barra };
839b072555Sjeremylt // Read q_data
849b072555Sjeremylt const CeedScalar w_det_J = q_data[i + Q * 0];
859b072555Sjeremylt // -- Grad-to-Grad q_data
86ed264d09SValeria Barra // ---- dXdx_j,k * dXdx_k,j
872b730f8bSJeremy L Thompson const CeedScalar dXdxdXdx_T[2][2] = {
882b730f8bSJeremy L Thompson {q_data[i + Q * 1], q_data[i + Q * 3]},
892b730f8bSJeremy L Thompson {q_data[i + Q * 3], q_data[i + Q * 2]}
90ed264d09SValeria Barra };
91ed264d09SValeria Barra
922b730f8bSJeremy L Thompson for (int k = 0; k < 3; k++) { // k = component
932b730f8bSJeremy L Thompson for (int j = 0; j < 2; j++) { // j = direction of vg
942b730f8bSJeremy L Thompson vJ[i + (k + j * 3) * Q] = w_det_J * (uJ[k][0] * dXdxdXdx_T[0][j] + uJ[k][1] * dXdxdXdx_T[1][j]);
952b730f8bSJeremy L Thompson }
962b730f8bSJeremy L Thompson }
97ed264d09SValeria Barra } // End of Quadrature Point Loop
98ed264d09SValeria Barra
99ed264d09SValeria Barra return 0;
100ed264d09SValeria Barra }
101ed264d09SValeria Barra // -----------------------------------------------------------------------------
102