// Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. // // SPDX-License-Identifier: BSD-2-Clause // // This file is part of CEED: http://github.com/ceed /// @file /// Operator for Navier-Stokes example using PETSc #ifndef newtonian_h #define newtonian_h #include #include #include "newtonian_state.h" #include "newtonian_types.h" #include "stabilization.h" #include "utils.h" // ***************************************************************************** // This QFunction sets a "still" initial condition for generic Newtonian IG problems // ***************************************************************************** CEED_QFUNCTION(ICsNewtonianIG)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { // Inputs const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; // Outputs CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; // Context const SetupContext context = (SetupContext)ctx; const CeedScalar theta0 = context->theta0; const CeedScalar P0 = context->P0; const CeedScalar cv = context->cv; const CeedScalar cp = context->cp; const CeedScalar *g = context->g; const CeedScalar Rd = cp - cv; // Quadrature Point Loop CeedPragmaSIMD for (CeedInt i=0; itheta0; const CeedScalar P0 = context->P0; // Quadrature Point Loop CeedPragmaSIMD for (CeedInt i=0; ig; const CeedScalar dt = context->dt; CeedPragmaSIMD // Quadrature Point Loop for (CeedInt i=0; ig; const CeedScalar dt = context->dt; CeedPragmaSIMD // Quadrature Point Loop for (CeedInt i=0; ig; CeedPragmaSIMD // Quadrature Point Loop for (CeedInt i=0; iijacobian_time_shift * dU[j] - dbody_force[j]); // -- Stabilization method: none (Galerkin), SU, or SUPG CeedScalar dstab[5][3], U_dot[5] = {0}; for (CeedInt j=0; j<5; j++) U_dot[j] = context->ijacobian_time_shift * dU[j]; Stabilization(context, s, Tau_d, grad_ds, U_dot, dbody_force, x_i, dstab); for (int j=0; j<5; j++) for (int k=0; k<3; k++) Grad_v[k][j][i] += wdetJ*(dstab[j][0] * dXdx[k][0] + dstab[j][1] * dXdx[k][1] + dstab[j][2] * dXdx[k][2]); } // End Quadrature Point Loop return 0; } CEED_QFUNCTION(IJacobian_Newtonian_Conserv)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return IJacobian_Newtonian(ctx, Q, in, out, StateFromU, StateFromU_fwd); } CEED_QFUNCTION(IJacobian_Newtonian_Prim)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return IJacobian_Newtonian(ctx, Q, in, out, StateFromY, StateFromY_fwd); } // ***************************************************************************** // Compute boundary integral (ie. for strongly set inflows) // ***************************************************************************** CEED_QFUNCTION_HELPER int BoundaryIntegral(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, StateFromQi_t StateFromQi, StateFromQi_fwd_t StateFromQi_fwd) { //*INDENT-OFF* const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*Grad_q)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], (*x)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[0], (*jac_data_sur)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA]) out[1]; //*INDENT-ON* const NewtonianIdealGasContext context = (NewtonianIdealGasContext) ctx; const bool is_implicit = context->is_implicit; CeedPragmaSIMD for(CeedInt i=0; iis_implicit; CeedPragmaSIMD // Quadrature Point Loop for (CeedInt i=0; iis_implicit; const CeedScalar P0 = context->P0; CeedPragmaSIMD // Quadrature Point Loop for (CeedInt i=0; iis_implicit; CeedPragmaSIMD // Quadrature Point Loop for (CeedInt i=0; iP0; ds.Y.pressure = 0.; State grad_ds[3]; for (CeedInt j=0; j<3; j++) { CeedScalar dx_i[3] = {0}, dqi_j[5]; for (CeedInt k=0; k<5; k++) dqi_j[k] = Grad_dq[0][k][i] * dXdx[0][j] + Grad_dq[1][k][i] * dXdx[1][j]; dx_i[j] = 1.; grad_ds[j] = StateFromQi_fwd(context, s, dqi_j, x_i, dx_i); } CeedScalar dstrain_rate[6], dkmstress[6], stress[3][3], dstress[3][3], dFe[3]; KMStrainRate(grad_ds, dstrain_rate); NewtonianStress(context, dstrain_rate, dkmstress); KMUnpack(dkmstress, dstress); KMUnpack(kmstress, stress); ViscousEnergyFlux_fwd(context, s.Y, ds.Y, grad_ds, stress, dstress, dFe); StateConservative dF_inviscid[3]; FluxInviscid_fwd(context, s, ds, dF_inviscid); CeedScalar dFlux[5]; FluxTotal_Boundary(dF_inviscid, dstress, dFe, norm, dFlux); for (int j=0; j<5; j++) v[j][i] = -wdetJb * dFlux[j]; } // End Quadrature Point Loop return 0; } CEED_QFUNCTION(PressureOutflow_Jacobian_Conserv)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return PressureOutflow_Jacobian(ctx, Q, in, out, StateFromU, StateFromU_fwd); } CEED_QFUNCTION(PressureOutflow_Jacobian_Prim)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { return PressureOutflow_Jacobian(ctx, Q, in, out, StateFromY, StateFromY_fwd); } #endif // newtonian_h