1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 33a8779fbSJames Wright 43a8779fbSJames Wright /// @file 53a8779fbSJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 63a8779fbSJames Wright 73a8779fbSJames Wright #include "../qfunctions/newtonian.h" 83a8779fbSJames Wright 9e419654dSJeremy L Thompson #include <ceed.h> 10e419654dSJeremy L Thompson #include <petscdm.h> 11e419654dSJeremy L Thompson 12149fb536SJames Wright #include <navierstokes.h> 136dd99bedSLeila Ghaffari 14e5a8cae0SJames Wright const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL}; 15e5a8cae0SJames Wright const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 16e5a8cae0SJames Wright 17cde3d787SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIGProperties gas); 18f0b65372SJed Brown 199b05e62eSJames Wright static PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) { 209b05e62eSJames Wright MPI_Comm comm = honee->comm; 219b05e62eSJames Wright Ceed ceed = honee->ceed; 22cde3d787SJames Wright NewtonianIdealGasContext newt_ctx; 239b05e62eSJames Wright 249b05e62eSJames Wright PetscFunctionBeginUser; 25cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newt_ctx)); 269b05e62eSJames Wright PetscCall(PetscPrintf(comm, 279b05e62eSJames Wright " Problem:\n" 289b05e62eSJames Wright " Problem Name : %s\n" 299b05e62eSJames Wright " Stabilization : %s\n", 30cde3d787SJames Wright app_ctx->problem_name, StabilizationTypes[newt_ctx->stabilization])); 31cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newt_ctx)); 329b05e62eSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 339b05e62eSJames Wright } 349b05e62eSJames Wright 3565dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 3665dee3d2SJames Wright // 3765dee3d2SJames Wright // Only used for SUPG stabilization 380c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(Honee honee, CeedOperator *op_mass) { 390c373b74SJames Wright Ceed ceed = honee->ceed; 4065dee3d2SJames Wright CeedInt num_comp_q, q_data_size; 4165dee3d2SJames Wright CeedQFunction qf_mass; 4201e19bfaSJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd; 4365dee3d2SJames Wright CeedBasis basis_q; 4465dee3d2SJames Wright CeedVector q_data; 45236a8ba6SJames Wright CeedQFunctionContext qfctx = NULL; 4665dee3d2SJames Wright PetscInt dim = 3; 4765dee3d2SJames Wright 4865dee3d2SJames Wright PetscFunctionBeginUser; 4965dee3d2SJames Wright { // Get restriction and basis from the RHS function 5065dee3d2SJames Wright CeedOperator *sub_ops; 5101e19bfaSJames Wright CeedOperatorField op_field; 5265dee3d2SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 5365dee3d2SJames Wright 54da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 5501e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 5601e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 5701e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 5801e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 5965dee3d2SJames Wright 60236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 6165dee3d2SJames Wright } 6265dee3d2SJames Wright 6365dee3d2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 6401e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 6565dee3d2SJames Wright 6665dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 6765dee3d2SJames Wright 68236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 6965dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 7065dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 7165dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 7265dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 7365dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 7465dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 7565dee3d2SJames Wright 7665dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 7771f2ed29SJames Wright PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Newtonian Stabilized")); 7865dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 790c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 8001e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 8165dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 8265dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 8365dee3d2SJames Wright 84fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 8501e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 86fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 87fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 88236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 8965dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 9065dee3d2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 9165dee3d2SJames Wright } 924c07ec22SJames Wright 9336038bbcSJames Wright /** 9436038bbcSJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 9536038bbcSJames Wright 960c373b74SJames Wright @param[in] honee `Honee` context 9736038bbcSJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 9836038bbcSJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 9936038bbcSJames Wright **/ 100e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 1010c373b74SJames Wright Ceed ceed = honee->ceed; 10236038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 103cf8f12c8SJames Wright PetscInt dim, num_comp_q; 104236a8ba6SJames Wright CeedQFunctionContext newtonian_qfctx = NULL; 105cf8f12c8SJames Wright CeedElemRestriction elem_restr_q; 106cf8f12c8SJames Wright CeedBasis basis_q; 10736038bbcSJames Wright 10836038bbcSJames Wright PetscFunctionBeginUser; 10936038bbcSJames Wright // -- Get Pre-requisite things 11036038bbcSJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 111cf8f12c8SJames Wright PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q)); 11236038bbcSJames Wright 1130880fbb6SJames Wright { // Get newtonian QF context 114b3b24828SJames Wright CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op; 11536038bbcSJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 11636038bbcSJames Wright 117da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(main_op, &sub_ops)); 118236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 11936038bbcSJames Wright } 120da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs)); 12136038bbcSJames Wright { // Add the volume integral CeedOperator 12236038bbcSJames Wright CeedQFunction qf_rhs_volume; 12336038bbcSJames Wright CeedOperator op_rhs_volume; 12436038bbcSJames Wright CeedVector q_data; 1250880fbb6SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 1260880fbb6SJames Wright CeedBasis basis_diff_flux = NULL; 12736038bbcSJames Wright CeedInt q_data_size; 12836038bbcSJames Wright 1290880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 130cf8f12c8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q)); 131cf8f12c8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q)); 132e3db12f8SJames Wright PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 133e3db12f8SJames Wright &elem_restr_qd, &q_data, &q_data_size)); 1340c373b74SJames Wright switch (honee->phys->state_var) { 13536038bbcSJames Wright case STATEVAR_PRIMITIVE: 13636038bbcSJames Wright PetscCallCeed(ceed, 13736038bbcSJames Wright CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume)); 13836038bbcSJames Wright break; 13936038bbcSJames Wright case STATEVAR_CONSERVATIVE: 14036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc, 14136038bbcSJames Wright &qf_rhs_volume)); 14236038bbcSJames Wright break; 14336038bbcSJames Wright case STATEVAR_ENTROPY: 14436038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc, 14536038bbcSJames Wright &qf_rhs_volume)); 14636038bbcSJames Wright break; 14736038bbcSJames Wright } 14836038bbcSJames Wright 149236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx)); 15036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP)); 15136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 15236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 15336038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 15436038bbcSJames Wright 15536038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 156cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 157cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 15836038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 15936038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 16036038bbcSJames Wright 161da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume)); 16236038bbcSJames Wright 16336038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 16436038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 1650880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 1660880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 167cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 168cf8f12c8SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 16936038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 17036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 17136038bbcSJames Wright } 17236038bbcSJames Wright 17336038bbcSJames Wright { // Add the boundary integral CeedOperator 17436038bbcSJames Wright CeedQFunction qf_rhs_boundary; 17536038bbcSJames Wright DMLabel face_sets_label; 17636038bbcSJames Wright PetscInt num_face_set_values, *face_set_values; 17736038bbcSJames Wright CeedInt q_data_size; 17836038bbcSJames Wright 17936038bbcSJames Wright // -- Build RHS operator 1800c373b74SJames Wright switch (honee->phys->state_var) { 18136038bbcSJames Wright case STATEVAR_PRIMITIVE: 18236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc, 18336038bbcSJames Wright &qf_rhs_boundary)); 18436038bbcSJames Wright break; 18536038bbcSJames Wright case STATEVAR_CONSERVATIVE: 18636038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc, 18736038bbcSJames Wright &qf_rhs_boundary)); 18836038bbcSJames Wright break; 18936038bbcSJames Wright case STATEVAR_ENTROPY: 19036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc, 19136038bbcSJames Wright &qf_rhs_boundary)); 19236038bbcSJames Wright break; 19336038bbcSJames Wright } 19436038bbcSJames Wright 1950c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 196236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx)); 19736038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP)); 19836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 19936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 20036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 20136038bbcSJames Wright 20236038bbcSJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 20336038bbcSJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 20436038bbcSJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 20536038bbcSJames Wright DMLabel face_orientation_label; 20636038bbcSJames Wright PetscInt num_orientations_values, *orientation_values; 20736038bbcSJames Wright 20836038bbcSJames Wright { 20936038bbcSJames Wright char *face_orientation_label_name; 21036038bbcSJames Wright 21136038bbcSJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 21236038bbcSJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 21336038bbcSJames Wright PetscCall(PetscFree(face_orientation_label_name)); 21436038bbcSJames Wright } 21536038bbcSJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 21636038bbcSJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 21736038bbcSJames Wright CeedOperator op_rhs_boundary; 21836038bbcSJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 21936038bbcSJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 22036038bbcSJames Wright CeedVector q_data; 22136038bbcSJames Wright CeedInt q_data_size; 22236038bbcSJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 22336038bbcSJames Wright 2240c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 2250c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 22636038bbcSJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 22736038bbcSJames Wright &elem_restr_diff_flux_boundary)); 2286b9fd993SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 2299eadbee4SJames Wright PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, 2309eadbee4SJames Wright &q_data_size)); 23136038bbcSJames Wright 23236038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 23336038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 23436038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 23536038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 23636038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 23736038bbcSJames Wright CEED_VECTOR_ACTIVE)); 23836038bbcSJames Wright 239da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_boundary)); 24036038bbcSJames Wright 24136038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 24236038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 24336038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 24436038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 24536038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 24636038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 24736038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 24836038bbcSJames Wright } 24936038bbcSJames Wright PetscCall(PetscFree(orientation_values)); 25036038bbcSJames Wright } 25136038bbcSJames Wright PetscCall(PetscFree(face_set_values)); 25236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 25336038bbcSJames Wright } 25436038bbcSJames Wright 255236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 25636038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 25736038bbcSJames Wright } 25836038bbcSJames Wright 25936038bbcSJames Wright /** 26036038bbcSJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 26136038bbcSJames Wright 2620c373b74SJames Wright @param[in] honee `Honee` context 26336038bbcSJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 26436038bbcSJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 26536038bbcSJames Wright **/ 266e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2670c373b74SJames Wright Ceed ceed = honee->ceed; 26836038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 269cf8f12c8SJames Wright CeedBasis basis_diff_flux, basis_q; 270cf8f12c8SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd, elem_restr_q; 27136038bbcSJames Wright CeedVector q_data; 272cf8f12c8SJames Wright CeedInt q_data_size; 273cf8f12c8SJames Wright PetscInt dim, num_comp_q; 274e3db12f8SJames Wright PetscInt height = 0, dm_field = 0; 27536038bbcSJames Wright CeedQFunction qf_rhs; 276236a8ba6SJames Wright CeedQFunctionContext newtonian_qfctx = NULL; 27736038bbcSJames Wright 27836038bbcSJames Wright PetscFunctionBeginUser; 27936038bbcSJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 280cf8f12c8SJames Wright PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q)); 28136038bbcSJames Wright 28240b78511SJames Wright { // Get newtonian QF context 283b3b24828SJames Wright CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op; 28436038bbcSJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 28536038bbcSJames Wright 286da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(main_op, &sub_ops)); 287236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 28836038bbcSJames Wright } 289cf8f12c8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q)); 290cf8f12c8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q)); 291e3db12f8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux)); 292e3db12f8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux)); 293e3db12f8SJames Wright PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 294e3db12f8SJames Wright &elem_restr_qd, &q_data, &q_data_size)); 29536038bbcSJames Wright 2960c373b74SJames Wright switch (honee->phys->state_var) { 29736038bbcSJames Wright case STATEVAR_PRIMITIVE: 29836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs)); 29936038bbcSJames Wright break; 30036038bbcSJames Wright case STATEVAR_CONSERVATIVE: 30136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs)); 30236038bbcSJames Wright break; 30336038bbcSJames Wright case STATEVAR_ENTROPY: 30436038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs)); 30536038bbcSJames Wright break; 30636038bbcSJames Wright } 30736038bbcSJames Wright 308236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx)); 30936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP)); 31036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 31136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 31236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 31336038bbcSJames Wright 31436038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 315cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 316cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 31736038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 31836038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 31936038bbcSJames Wright 32036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 321236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 32236038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 32336038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 32436038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 32536038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 326cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 327cf8f12c8SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 32836038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 32936038bbcSJames Wright } 33036038bbcSJames Wright 331d6cac220SJames Wright static PetscErrorCode BoundaryIntegralBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 332d6cac220SJames Wright HoneeBCStruct honee_bc; 333d6cac220SJames Wright 334d6cac220SJames Wright PetscFunctionBeginUser; 335d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 336d6cac220SJames Wright Honee honee = honee_bc->honee; 337d6cac220SJames Wright 338d6cac220SJames Wright switch (honee->phys->state_var) { 339d6cac220SJames Wright case STATEVAR_CONSERVATIVE: 340d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Conserv, BoundaryIntegral_Conserv_loc, honee_bc->qfctx, qf)); 341d6cac220SJames Wright break; 342d6cac220SJames Wright case STATEVAR_PRIMITIVE: 343d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Prim, BoundaryIntegral_Prim_loc, honee_bc->qfctx, qf)); 344d6cac220SJames Wright break; 345d6cac220SJames Wright case STATEVAR_ENTROPY: 346d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Entropy, BoundaryIntegral_Entropy_loc, honee_bc->qfctx, qf)); 347d6cac220SJames Wright break; 348d6cac220SJames Wright } 349d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 350d6cac220SJames Wright } 351d6cac220SJames Wright 352d6cac220SJames Wright static PetscErrorCode BoundaryIntegralBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) { 353d6cac220SJames Wright HoneeBCStruct honee_bc; 354d6cac220SJames Wright 355d6cac220SJames Wright PetscFunctionBeginUser; 356d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 357d6cac220SJames Wright Honee honee = honee_bc->honee; 358d6cac220SJames Wright 359d6cac220SJames Wright switch (honee->phys->state_var) { 360d6cac220SJames Wright case STATEVAR_CONSERVATIVE: 361d6cac220SJames Wright PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Conserv, BoundaryIntegral_Jacobian_Conserv_loc, honee_bc->qfctx, qf)); 362d6cac220SJames Wright break; 363d6cac220SJames Wright case STATEVAR_PRIMITIVE: 364d6cac220SJames Wright PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Prim, BoundaryIntegral_Jacobian_Prim_loc, honee_bc->qfctx, qf)); 365d6cac220SJames Wright break; 366d6cac220SJames Wright case STATEVAR_ENTROPY: 367d6cac220SJames Wright PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Entropy, BoundaryIntegral_Jacobian_Entropy_loc, honee_bc->qfctx, qf)); 368d6cac220SJames Wright break; 369d6cac220SJames Wright } 370d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 371d6cac220SJames Wright } 372d6cac220SJames Wright 373d3c60affSJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx) { 374b7f03f12SJed Brown SetupContext setup_context; 3750c373b74SJames Wright Honee honee = *(Honee *)ctx; 3760c373b74SJames Wright CeedInt degree = honee->app_ctx->degree; 3773a8779fbSJames Wright StabilizationType stab; 3783636f6a4SJames Wright StateVariable state_var; 3790c373b74SJames Wright MPI_Comm comm = honee->comm; 3800c373b74SJames Wright Ceed ceed = honee->ceed; 3813a8779fbSJames Wright PetscBool implicit; 3820d8cd818SJames Wright PetscBool unit_tests; 38315a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 3843a8779fbSJames Wright 38515a3537eSJed Brown PetscFunctionBeginUser; 386bff204bcSJames Wright // Option Defaults 387f5dc303cSJames Wright const CeedScalar Cv_func[3] = {36, 60, 128}; 388bff204bcSJames Wright StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 389bff204bcSJames Wright CeedScalar idl_decay_time = -1; 390bff204bcSJames Wright PetscCall(PetscNew(&newtonian_ig_ctx)); 391f5dc303cSJames Wright *newtonian_ig_ctx = (struct NewtonianIdealGasContext_){ 392cde3d787SJames Wright .gas = 393cde3d787SJames Wright { 394cde3d787SJames Wright .cv = 717., 395cde3d787SJames Wright .cp = 1004., 396cde3d787SJames Wright .lambda = -2. / 3., 397cde3d787SJames Wright .mu = 1.8e-5, 398cde3d787SJames Wright .k = 0.02638, 399cde3d787SJames Wright }, 400cde3d787SJames Wright .tau_coeffs = 401cde3d787SJames Wright { 402cde3d787SJames Wright .Ctau_t = 1.0, 403bff204bcSJames Wright .Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1], 404bff204bcSJames Wright .Ctau_C = 0.25 / degree, 405bff204bcSJames Wright .Ctau_M = 0.25 / degree, 406bff204bcSJames Wright .Ctau_E = 0.125, 407cde3d787SJames Wright }, 408cde3d787SJames Wright .g = {0, 0, 0}, // m/s^2 409bff204bcSJames Wright .idl_start = 0, 410bff204bcSJames Wright .idl_length = 0, 411bff204bcSJames Wright .idl_pressure = reference.pressure, 412bff204bcSJames Wright .idl_enable = PETSC_FALSE, 413bff204bcSJames Wright }; 414bff204bcSJames Wright 4158d8b1f0cSJames Wright PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 4168d8b1f0cSJames Wright PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 4178d8b1f0cSJames Wright (PetscEnum *)&state_var, NULL)); 4183a8779fbSJames Wright 4198d8b1f0cSJames Wright // Newtonian fluid properties 420cde3d787SJames Wright PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, newtonian_ig_ctx->gas.cv, &newtonian_ig_ctx->gas.cv, NULL)); 421cde3d787SJames Wright PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, newtonian_ig_ctx->gas.cp, &newtonian_ig_ctx->gas.cp, NULL)); 422cde3d787SJames Wright PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, newtonian_ig_ctx->gas.lambda, 423cde3d787SJames Wright &newtonian_ig_ctx->gas.lambda, NULL)); 424cde3d787SJames Wright PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, newtonian_ig_ctx->gas.mu, &newtonian_ig_ctx->gas.mu, NULL)); 425cde3d787SJames Wright PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, newtonian_ig_ctx->gas.k, &newtonian_ig_ctx->gas.k, NULL)); 4268d8b1f0cSJames Wright 4278d8b1f0cSJames Wright PetscInt dim = 3; 428bff204bcSJames Wright PetscBool given_option = PETSC_FALSE; 4298d8b1f0cSJames Wright PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 430bff204bcSJames Wright PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, newtonian_ig_ctx->g, &dim, &given_option)); 4318d8b1f0cSJames Wright if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 4328d8b1f0cSJames Wright 4338d8b1f0cSJames Wright // Stabilization parameters 4348d8b1f0cSJames Wright PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 435cde3d787SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_t, 436cde3d787SJames Wright &newtonian_ig_ctx->tau_coeffs.Ctau_t, NULL)); 437cde3d787SJames Wright PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_v, 438cde3d787SJames Wright &newtonian_ig_ctx->tau_coeffs.Ctau_v, NULL)); 439cde3d787SJames Wright PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_C, 440cde3d787SJames Wright &newtonian_ig_ctx->tau_coeffs.Ctau_C, NULL)); 441cde3d787SJames Wright PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_M, 442cde3d787SJames Wright &newtonian_ig_ctx->tau_coeffs.Ctau_M, NULL)); 443cde3d787SJames Wright PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_E, 444cde3d787SJames Wright &newtonian_ig_ctx->tau_coeffs.Ctau_E, NULL)); 4458d8b1f0cSJames Wright 4468d8b1f0cSJames Wright dim = 3; 4478d8b1f0cSJames Wright PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 4488d8b1f0cSJames Wright PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 4498d8b1f0cSJames Wright PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 4508d8b1f0cSJames Wright 4518d8b1f0cSJames Wright PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 4528d8b1f0cSJames Wright PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 4538d8b1f0cSJames Wright PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 4548d8b1f0cSJames Wright "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 4558d8b1f0cSJames Wright 4568d8b1f0cSJames Wright // IDL Settings 457bff204bcSJames Wright { 458bff204bcSJames Wright PetscBool idl_enable = (PetscBool)newtonian_ig_ctx->idl_enable; // Need PetscBool variable to read in from PetscOptionsScalar() 4598d8b1f0cSJames Wright PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 4608d8b1f0cSJames Wright NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 4618d8b1f0cSJames Wright PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 4628d8b1f0cSJames Wright if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 463bff204bcSJames Wright newtonian_ig_ctx->idl_enable = idl_enable; 464bff204bcSJames Wright 4659eadbee4SJames Wright PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, newtonian_ig_ctx->idl_start, &newtonian_ig_ctx->idl_start, 4669eadbee4SJames Wright NULL)); 467bff204bcSJames Wright PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, newtonian_ig_ctx->idl_length, 468bff204bcSJames Wright &newtonian_ig_ctx->idl_length, NULL)); 469bff204bcSJames Wright newtonian_ig_ctx->idl_pressure = reference.pressure; 470bff204bcSJames Wright PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, 471bff204bcSJames Wright newtonian_ig_ctx->idl_pressure, &newtonian_ig_ctx->idl_pressure, NULL)); 472bff204bcSJames Wright } 4738d8b1f0cSJames Wright PetscOptionsEnd(); 4748d8b1f0cSJames Wright 475f5dc303cSJames Wright // ------------------------------------------------------ 476f5dc303cSJames Wright // Set up the QFunction context 477f5dc303cSJames Wright // ------------------------------------------------------ 478f5dc303cSJames Wright // -- Scale variables to desired units 479f5dc303cSJames Wright Units units = honee->units; 480cde3d787SJames Wright newtonian_ig_ctx->gas.cv *= units->J_per_kg_K; 481cde3d787SJames Wright newtonian_ig_ctx->gas.cp *= units->J_per_kg_K; 482cde3d787SJames Wright newtonian_ig_ctx->gas.mu *= units->Pascal * units->second; 483cde3d787SJames Wright newtonian_ig_ctx->gas.k *= units->W_per_m_K; 484f5dc303cSJames Wright for (PetscInt i = 0; i < 3; i++) newtonian_ig_ctx->g[i] *= units->m_per_squared_s; 485f5dc303cSJames Wright reference.pressure *= units->Pascal; 486f5dc303cSJames Wright for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= units->meter / units->second; 487f5dc303cSJames Wright reference.temperature *= units->Kelvin; 488f5dc303cSJames Wright 489f5dc303cSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 490f5dc303cSJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 491*7e3656bdSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = (domain_max[i] - domain_min[i]); 492f5dc303cSJames Wright 493f5dc303cSJames Wright // -- Solver Settings 494f5dc303cSJames Wright honee->phys->implicit = implicit; 495f5dc303cSJames Wright honee->phys->state_var = state_var; 496f5dc303cSJames Wright 497f5dc303cSJames Wright // -- QFunction Context 498f5dc303cSJames Wright newtonian_ig_ctx->stabilization = stab; 499f5dc303cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 500f5dc303cSJames Wright newtonian_ig_ctx->state_var = state_var; 501f5dc303cSJames Wright newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * units->second); 502f5dc303cSJames Wright newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 503f5dc303cSJames Wright 504f5dc303cSJames Wright // -- Setup Context 505f5dc303cSJames Wright PetscCall(PetscNew(&setup_context)); 506f5dc303cSJames Wright *setup_context = (struct SetupContext_){ 507f5dc303cSJames Wright .reference = reference, 508cde3d787SJames Wright .newt_ctx = *newtonian_ig_ctx, 509f5dc303cSJames Wright .lx = domain_size[0], 510f5dc303cSJames Wright .ly = domain_size[1], 511f5dc303cSJames Wright .lz = domain_size[2], 512f5dc303cSJames Wright .time = 0, 513f5dc303cSJames Wright }; 514f5dc303cSJames Wright 515f5dc303cSJames Wright CeedQFunctionContext ics_qfctx, newtonian_ig_qfctx; 516f5dc303cSJames Wright 517f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); 518f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 519f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 520f5dc303cSJames Wright PetscCallCeed(ceed, 521f5dc303cSJames Wright CeedQFunctionContextRegisterDouble(ics_qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation")); 522f5dc303cSJames Wright 523f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx)); 524f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 525f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 526f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 527f5dc303cSJames Wright "Size of timestep, delta t")); 528f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift", 529f5dc303cSJames Wright offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 530f5dc303cSJames Wright "Shift for mass matrix in IJacobian")); 531f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 532f5dc303cSJames Wright "Current solution time")); 533f5dc303cSJames Wright 5348d8b1f0cSJames Wright // Set problem information 5351abc2837SJames Wright problem->num_comps_jac_data = 14; 536bff204bcSJames Wright if (newtonian_ig_ctx->idl_enable) problem->num_comps_jac_data += 1; 53758ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 538cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 539f27c2204SJames Wright problem->num_components = 5; 540f27c2204SJames Wright PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 541f5dc303cSJames Wright static const char *const conserv_component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 542f5dc303cSJames Wright static const char *const prim_component_names[] = {"Pressure", "VelocityX", "VelocityY", "VelocityZ", "Temperature"}; 5438d8b1f0cSJames Wright static const char *const entropy_component_names[] = {"EntropyDensity", "EntropyMomentumX", "EntropyMomentumY", "EntropyMomentumZ", 5448d8b1f0cSJames Wright "EntropyTotalEnergy"}; 5458d8b1f0cSJames Wright 5468d8b1f0cSJames Wright switch (state_var) { 547f27c2204SJames Wright case STATEVAR_CONSERVATIVE: 548f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Conserv, .qf_loc = ICsNewtonianIG_Conserv_loc}; 549f5dc303cSJames Wright problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHSFunction_Newtonian, .qf_loc = RHSFunction_Newtonian_loc}; 550f5dc303cSJames Wright problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Conserv, .qf_loc = IFunction_Newtonian_Conserv_loc}; 551f5dc303cSJames Wright problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Conserv, .qf_loc = IJacobian_Newtonian_Conserv_loc}; 552f5dc303cSJames Wright for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(conserv_component_names[i], &problem->component_names[i])); 553f27c2204SJames Wright break; 554f27c2204SJames Wright case STATEVAR_PRIMITIVE: 555f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Prim, .qf_loc = ICsNewtonianIG_Prim_loc}; 556f5dc303cSJames Wright problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Prim, .qf_loc = IFunction_Newtonian_Prim_loc}; 557f5dc303cSJames Wright problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Prim, .qf_loc = IJacobian_Newtonian_Prim_loc}; 558f5dc303cSJames Wright for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(prim_component_names[i], &problem->component_names[i])); 559f27c2204SJames Wright break; 560f27c2204SJames Wright case STATEVAR_ENTROPY: 561f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Entropy, .qf_loc = ICsNewtonianIG_Entropy_loc}; 562f5dc303cSJames Wright problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Entropy, .qf_loc = IFunction_Newtonian_Entropy_loc}; 563f5dc303cSJames Wright problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Entropy, .qf_loc = IJacobian_Newtonian_Entropy_loc}; 564f27c2204SJames Wright for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(entropy_component_names[i], &problem->component_names[i])); 565f27c2204SJames Wright break; 566f27c2204SJames Wright } 567f5dc303cSJames Wright // All QFunctions get the same QFunctionContext regardless of state variable 568f5dc303cSJames Wright problem->ics.qfctx = ics_qfctx; 569f5dc303cSJames Wright problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; 570f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx)); 571f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx)); 572f27c2204SJames Wright 5738d8b1f0cSJames Wright if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 5748d8b1f0cSJames Wright 575d9e69621SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 4, &honee->diff_flux_proj)); 5760c373b74SJames Wright if (honee->diff_flux_proj) { 5770c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 57836038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 579ae2e5212SJames Wright PetscSection section; 58036038bbcSJames Wright 58136038bbcSJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_NS; 58236038bbcSJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS; 583ae2e5212SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 5840c373b74SJames Wright switch (honee->diff_flux_proj->method) { 58536038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 58636038bbcSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 58736038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX")); 58836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY")); 58936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ")); 59036038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy")); 59136038bbcSJames Wright } break; 59236038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 59336038bbcSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 59436038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX")); 59536038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY")); 59636038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ")); 59736038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX")); 59836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY")); 59936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ")); 60036038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX")); 60136038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY")); 60236038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ")); 60336038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX")); 60436038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY")); 60536038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ")); 60636038bbcSJames Wright } break; 60736038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 6080c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 6090c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 61036038bbcSJames Wright break; 61136038bbcSJames Wright } 61236038bbcSJames Wright } 61336038bbcSJames Wright 6145e79d562SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 6155e79d562SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 6165e79d562SJames Wright const char *name; 6175e79d562SJames Wright 6185e79d562SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 6195e79d562SJames Wright if (!strcmp(name, "slip")) { 6205e79d562SJames Wright PetscCall(SlipBCSetup(bc_def, problem, dm, ctx, newtonian_ig_qfctx)); 621d4e423e7SJames Wright } else if (!strcmp(name, "freestream")) { 622d4e423e7SJames Wright PetscCall(FreestreamBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference)); 623f978755dSJames Wright } else if (!strcmp(name, "outflow")) { 624f978755dSJames Wright PetscCall(OutflowBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference)); 625d6cac220SJames Wright } else if (!strcmp(name, "inflow")) { 626d6cac220SJames Wright HoneeBCStruct honee_bc; 627d6cac220SJames Wright 628d6cac220SJames Wright PetscCall(PetscNew(&honee_bc)); 629d6cac220SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &honee_bc->qfctx)); 630d6cac220SJames Wright honee_bc->honee = honee; 6311abc2837SJames Wright honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0; 632d6cac220SJames Wright PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 633d6cac220SJames Wright 634d6cac220SJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, BoundaryIntegralBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 635d6cac220SJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, BoundaryIntegralBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp)); 6365e79d562SJames Wright } 6375e79d562SJames Wright } 6389ed3d70dSJames Wright 639cde3d787SJames Wright if (unit_tests) PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx->gas)); 640d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6413a8779fbSJames Wright } 6423a8779fbSJames Wright 6439b05e62eSJames Wright //------------------------------------ 6449b05e62eSJames Wright // Unit test functions 6459b05e62eSJames Wright //------------------------------------ 64648417fb7SJames Wright 64748417fb7SJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name, 64848417fb7SJames Wright PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) { 64948417fb7SJames Wright CeedScalar relative_error[5]; // relative error 65048417fb7SJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 65148417fb7SJames Wright 65248417fb7SJames Wright PetscFunctionBeginUser; 65348417fb7SJames Wright relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1); 65448417fb7SJames Wright relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1); 65548417fb7SJames Wright 65648417fb7SJames Wright CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3])); 65748417fb7SJames Wright CeedScalar u_divisor = u_magnitude > divisor_threshold ? u_magnitude : 1; 65848417fb7SJames Wright for (int i = 1; i < 4; i++) { 65948417fb7SJames Wright relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor; 66048417fb7SJames Wright } 66148417fb7SJames Wright 66248417fb7SJames Wright if (fabs(relative_error[0]) >= rtol_0) { 66348417fb7SJames Wright printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]); 66448417fb7SJames Wright } 66548417fb7SJames Wright for (int i = 1; i < 4; i++) { 66648417fb7SJames Wright if (fabs(relative_error[i]) >= rtol_u) { 66748417fb7SJames Wright printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]); 66848417fb7SJames Wright } 66948417fb7SJames Wright } 67048417fb7SJames Wright if (fabs(relative_error[4]) >= rtol_4) { 67148417fb7SJames Wright printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]); 67248417fb7SJames Wright } 67348417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 67448417fb7SJames Wright } 67548417fb7SJames Wright 67648417fb7SJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test 677cde3d787SJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIGProperties gas, const CeedScalar A0[5], 67848417fb7SJames Wright CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 67948417fb7SJames Wright CeedScalar B0[5], A0_test[5]; 68048417fb7SJames Wright char buf[128]; 68148417fb7SJames Wright const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 68248417fb7SJames Wright 68348417fb7SJames Wright PetscFunctionBeginUser; 68448417fb7SJames Wright const char *A_initial = StateVariables_Initial[state_var_A]; 68548417fb7SJames Wright const char *B_initial = StateVariables_Initial[state_var_B]; 68648417fb7SJames Wright 68748417fb7SJames Wright State state_A0 = StateFromQ(gas, A0, state_var_A); 68848417fb7SJames Wright StateToQ(gas, state_A0, B0, state_var_B); 68948417fb7SJames Wright State state_B0 = StateFromQ(gas, B0, state_var_B); 69048417fb7SJames Wright StateToQ(gas, state_B0, A0_test, state_var_A); 69148417fb7SJames Wright 69248417fb7SJames Wright snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial); 69348417fb7SJames Wright PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4)); 69448417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 69548417fb7SJames Wright } 69648417fb7SJames Wright 69748417fb7SJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation 698cde3d787SJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIGProperties gas, const CeedScalar A0[5], 69948417fb7SJames Wright CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 70048417fb7SJames Wright CeedScalar eps = 4e-7; // Finite difference step 70148417fb7SJames Wright char buf[128]; 70248417fb7SJames Wright const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 70348417fb7SJames Wright 70448417fb7SJames Wright PetscFunctionBeginUser; 70548417fb7SJames Wright const char *A_initial = StateVariables_Initial[state_var_A]; 70648417fb7SJames Wright const char *B_initial = StateVariables_Initial[state_var_B]; 70748417fb7SJames Wright State state_0 = StateFromQ(gas, A0, state_var_A); 70848417fb7SJames Wright 70948417fb7SJames Wright for (int i = 0; i < 5; i++) { 71048417fb7SJames Wright CeedScalar dB[5] = {0.}, dB_fd[5] = {0.}; 71148417fb7SJames Wright { // Calculate dB using State functions 71248417fb7SJames Wright CeedScalar dA[5] = {0}; 71348417fb7SJames Wright 71448417fb7SJames Wright dA[i] = A0[i]; 71548417fb7SJames Wright State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A); 71648417fb7SJames Wright StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B); 71748417fb7SJames Wright } 71848417fb7SJames Wright 71948417fb7SJames Wright { // Calculate dB_fd via finite difference approximation 72048417fb7SJames Wright CeedScalar A1[5], B0[5], B1[5]; 72148417fb7SJames Wright 72248417fb7SJames Wright for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j]; 72348417fb7SJames Wright State state_1 = StateFromQ(gas, A1, state_var_A); 72448417fb7SJames Wright StateToQ(gas, state_0, B0, state_var_B); 72548417fb7SJames Wright StateToQ(gas, state_1, B1, state_var_B); 72648417fb7SJames Wright for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps; 72748417fb7SJames Wright } 72848417fb7SJames Wright 72948417fb7SJames Wright snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial); 73048417fb7SJames Wright PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4)); 73148417fb7SJames Wright } 73248417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 73348417fb7SJames Wright } 73448417fb7SJames Wright 73548417fb7SJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*` 736cde3d787SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIGProperties gas) { 73748417fb7SJames Wright Units units = honee->units; 73848417fb7SJames Wright const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin; 73948417fb7SJames Wright CeedScalar rtol; 74048417fb7SJames Wright 74148417fb7SJames Wright PetscFunctionBeginUser; 74248417fb7SJames Wright const CeedScalar T = 200 * K; 74348417fb7SJames Wright const CeedScalar rho = 1.2 * kg / Cube(m); 744cde3d787SJames Wright const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas.cv * T; 74548417fb7SJames Wright const CeedScalar u_base = 40 * m / sec; 74648417fb7SJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 74748417fb7SJames Wright const CeedScalar e_kinetic = 0.5 * Dot3(u, u); 748cde3d787SJames Wright const CeedScalar e_internal = gas.cv * T; 74948417fb7SJames Wright const CeedScalar e_total = e_kinetic + e_internal; 75048417fb7SJames Wright const CeedScalar gamma = HeatCapacityRatio(gas); 75148417fb7SJames Wright const CeedScalar entropy = log(P) - gamma * log(rho); 75248417fb7SJames Wright const CeedScalar rho_div_p = rho / P; 75348417fb7SJames Wright const CeedScalar Y0[5] = {P, u[0], u[1], u[2], T}; 75448417fb7SJames Wright const CeedScalar U0[5] = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total}; 75548417fb7SJames Wright const CeedScalar V0[5] = {(gamma - entropy) / (gamma - 1) - rho_div_p * (e_kinetic), rho_div_p * u[0], rho_div_p * u[1], rho_div_p * u[2], 75648417fb7SJames Wright -rho_div_p}; 75748417fb7SJames Wright 75848417fb7SJames Wright rtol = 20 * CEED_EPSILON; 75948417fb7SJames Wright PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 76048417fb7SJames Wright PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 76148417fb7SJames Wright PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 76248417fb7SJames Wright PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol)); 76348417fb7SJames Wright PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol)); 76448417fb7SJames Wright PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol)); 76548417fb7SJames Wright 76648417fb7SJames Wright rtol = 5e-6; 76748417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 76848417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 76948417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 77048417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol)); 77148417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol)); 77248417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol)); 77348417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 77448417fb7SJames Wright } 775