1*3c11f1fcSJeremy L Thompson // Copyright (c) 2017-2025, Lawrence Livermore National Security, LLC and other CEED contributors. 2*3c11f1fcSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3*3c11f1fcSJeremy L Thompson // 4*3c11f1fcSJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5*3c11f1fcSJeremy L Thompson // 6*3c11f1fcSJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7*3c11f1fcSJeremy L Thompson 8*3c11f1fcSJeremy L Thompson /// @file 9*3c11f1fcSJeremy L Thompson /// libCEED QFunctions for diffusion operator example using PETSc 10*3c11f1fcSJeremy L Thompson 11*3c11f1fcSJeremy L Thompson #include <ceed/types.h> 12*3c11f1fcSJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS 13*3c11f1fcSJeremy L Thompson #include <math.h> 14*3c11f1fcSJeremy L Thompson #endif 15*3c11f1fcSJeremy L Thompson 16*3c11f1fcSJeremy L Thompson // ----------------------------------------------------------------------------- 17*3c11f1fcSJeremy L Thompson // This QFunction sets up the rhs and true solution for the problem 18*3c11f1fcSJeremy L Thompson // ----------------------------------------------------------------------------- 19*3c11f1fcSJeremy L Thompson CEED_QFUNCTION(SetupMassDiffRhs3)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 20*3c11f1fcSJeremy L Thompson #ifndef M_PI 21*3c11f1fcSJeremy L Thompson #define M_PI 3.14159265358979323846 22*3c11f1fcSJeremy L Thompson #endif 23*3c11f1fcSJeremy L Thompson const CeedScalar *x = in[0], *w = in[1]; 24*3c11f1fcSJeremy L Thompson CeedScalar *true_soln = out[0], *rhs = out[1]; 25*3c11f1fcSJeremy L Thompson 26*3c11f1fcSJeremy L Thompson // Quadrature Point Loop 27*3c11f1fcSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 28*3c11f1fcSJeremy L Thompson const CeedScalar c[3] = {0, 1., 2.}; 29*3c11f1fcSJeremy L Thompson const CeedScalar k[3] = {1., 2., 3.}; 30*3c11f1fcSJeremy L Thompson 31*3c11f1fcSJeremy L Thompson // Component 1 32*3c11f1fcSJeremy L Thompson true_soln[i + 0 * Q] = 33*3c11f1fcSJeremy L Thompson sin(M_PI * (c[0] + k[0] * x[i + Q * 0])) * sin(M_PI * (c[1] + k[1] * x[i + Q * 1])) * sin(M_PI * (c[2] + k[2] * x[i + Q * 2])); 34*3c11f1fcSJeremy L Thompson // Component 2 35*3c11f1fcSJeremy L Thompson true_soln[i + 1 * Q] = 2 * true_soln[i + 0 * Q]; 36*3c11f1fcSJeremy L Thompson // Component 3 37*3c11f1fcSJeremy L Thompson true_soln[i + 2 * Q] = 3 * true_soln[i + 0 * Q]; 38*3c11f1fcSJeremy L Thompson 39*3c11f1fcSJeremy L Thompson // Component 1 40*3c11f1fcSJeremy L Thompson rhs[i + 0 * Q] = w[i + Q * 0] * (M_PI * M_PI * (k[0] * k[0] + k[1] * k[1] + k[2] * k[2]) + 1.0) * true_soln[i + 0 * Q]; 41*3c11f1fcSJeremy L Thompson // Component 2 42*3c11f1fcSJeremy L Thompson rhs[i + 1 * Q] = 2 * rhs[i + 0 * Q]; 43*3c11f1fcSJeremy L Thompson // Component 3 44*3c11f1fcSJeremy L Thompson rhs[i + 2 * Q] = 3 * rhs[i + 0 * Q]; 45*3c11f1fcSJeremy L Thompson } // End of Quadrature Point Loop 46*3c11f1fcSJeremy L Thompson return 0; 47*3c11f1fcSJeremy L Thompson } 48*3c11f1fcSJeremy L Thompson 49*3c11f1fcSJeremy L Thompson // ----------------------------------------------------------------------------- 50*3c11f1fcSJeremy L Thompson // This QFunction applies the mass + diffusion operator for a vector field of 3 components. 51*3c11f1fcSJeremy L Thompson // 52*3c11f1fcSJeremy L Thompson // Inputs: 53*3c11f1fcSJeremy L Thompson // u - Input vector at quadrature points 54*3c11f1fcSJeremy L Thompson // ug - Input vector Jacobian at quadrature points 55*3c11f1fcSJeremy L Thompson // q_data - Geometric factors 56*3c11f1fcSJeremy L Thompson // 57*3c11f1fcSJeremy L Thompson // Output: 58*3c11f1fcSJeremy L Thompson // v - Output vector (test functions) at quadrature points 59*3c11f1fcSJeremy L Thompson // vJ - Output vector (test functions) Jacobian at quadrature points 60*3c11f1fcSJeremy L Thompson // ----------------------------------------------------------------------------- 61*3c11f1fcSJeremy L Thompson CEED_QFUNCTION(MassDiff3)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 62*3c11f1fcSJeremy L Thompson const CeedScalar *u = in[0], *ug = in[1], *q_data = in[2]; 63*3c11f1fcSJeremy L Thompson CeedScalar *v = out[0], *vg = out[1]; 64*3c11f1fcSJeremy L Thompson 65*3c11f1fcSJeremy L Thompson // Quadrature Point Loop 66*3c11f1fcSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 67*3c11f1fcSJeremy L Thompson // Read spatial derivatives of u components 68*3c11f1fcSJeremy L Thompson const CeedScalar uJ[3][3] = { 69*3c11f1fcSJeremy L Thompson {ug[i + (0 + 0 * 3) * Q], ug[i + (0 + 1 * 3) * Q], ug[i + (0 + 2 * 3) * Q]}, 70*3c11f1fcSJeremy L Thompson {ug[i + (1 + 0 * 3) * Q], ug[i + (1 + 1 * 3) * Q], ug[i + (1 + 2 * 3) * Q]}, 71*3c11f1fcSJeremy L Thompson {ug[i + (2 + 0 * 3) * Q], ug[i + (2 + 1 * 3) * Q], ug[i + (2 + 2 * 3) * Q]} 72*3c11f1fcSJeremy L Thompson }; 73*3c11f1fcSJeremy L Thompson // Read q_data (dXdxdXdx_T symmetric matrix) 74*3c11f1fcSJeremy L Thompson const CeedScalar dXdxdXdx_T[3][3] = { 75*3c11f1fcSJeremy L Thompson {q_data[i + 1 * Q], q_data[i + 2 * Q], q_data[i + 3 * Q]}, 76*3c11f1fcSJeremy L Thompson {q_data[i + 2 * Q], q_data[i + 4 * Q], q_data[i + 5 * Q]}, 77*3c11f1fcSJeremy L Thompson {q_data[i + 3 * Q], q_data[i + 5 * Q], q_data[i + 6 * Q]} 78*3c11f1fcSJeremy L Thompson }; 79*3c11f1fcSJeremy L Thompson 80*3c11f1fcSJeremy L Thompson for (int k = 0; k < 3; k++) { // k = component 81*3c11f1fcSJeremy L Thompson // Mass 82*3c11f1fcSJeremy L Thompson v[i + k * Q] = q_data[i + 0 * Q] * u[i + k * Q]; 83*3c11f1fcSJeremy L Thompson // Diff 84*3c11f1fcSJeremy L Thompson for (int j = 0; j < 3; j++) { // j = direction of vg 85*3c11f1fcSJeremy L Thompson vg[i + (k + j * 3) * Q] = (uJ[k][0] * dXdxdXdx_T[0][j] + uJ[k][1] * dXdxdXdx_T[1][j] + uJ[k][2] * dXdxdXdx_T[2][j]); 86*3c11f1fcSJeremy L Thompson } 87*3c11f1fcSJeremy L Thompson } 88*3c11f1fcSJeremy L Thompson } // End of Quadrature Point Loop 89*3c11f1fcSJeremy L Thompson 90*3c11f1fcSJeremy L Thompson return 0; 91*3c11f1fcSJeremy L Thompson } 92*3c11f1fcSJeremy L Thompson // ----------------------------------------------------------------------------- 93