1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utility functions for setting up problems using the Newtonian Qfunction 6 7 #include "../qfunctions/newtonian.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 12 #include <navierstokes.h> 13 14 const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL}; 15 const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 16 17 static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIGProperties gas); 18 19 static PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) { 20 MPI_Comm comm = honee->comm; 21 Ceed ceed = honee->ceed; 22 NewtonianIdealGasContext newt_ctx; 23 24 PetscFunctionBeginUser; 25 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newt_ctx)); 26 PetscCall(PetscPrintf(comm, 27 " Problem:\n" 28 " Problem Name : %s\n" 29 " Stabilization : %s\n", 30 app_ctx->problem_name, StabilizationTypes[newt_ctx->stabilization])); 31 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newt_ctx)); 32 PetscFunctionReturn(PETSC_SUCCESS); 33 } 34 35 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 36 // 37 // Only used for SUPG stabilization 38 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(Honee honee, CeedOperator *op_mass) { 39 Ceed ceed = honee->ceed; 40 CeedInt num_comp_q, q_data_size; 41 CeedQFunction qf_mass; 42 CeedElemRestriction elem_restr_q, elem_restr_qd; 43 CeedBasis basis_q; 44 CeedVector q_data; 45 CeedQFunctionContext qfctx = NULL; 46 PetscInt dim = 3; 47 48 PetscFunctionBeginUser; 49 { // Get restriction and basis from the RHS function 50 CeedOperator *sub_ops; 51 CeedOperatorField op_field; 52 PetscInt sub_op_index = 0; // will be 0 for the volume op 53 54 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 55 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 56 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 57 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 58 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 59 60 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 61 } 62 63 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 64 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 65 66 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 67 68 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 69 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 70 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 71 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 72 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 73 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 74 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 75 76 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 77 PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Newtonian Stabilized")); 78 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 79 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 80 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 81 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 82 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 83 84 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 85 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 86 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 87 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 88 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 89 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 90 PetscFunctionReturn(PETSC_SUCCESS); 91 } 92 93 /** 94 @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 95 96 @param[in] honee `Honee` context 97 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 98 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 99 **/ 100 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 101 Ceed ceed = honee->ceed; 102 NodalProjectionData projection = diff_flux_proj->projection; 103 PetscInt dim, num_comp_q; 104 CeedQFunctionContext newtonian_qfctx = NULL; 105 CeedElemRestriction elem_restr_q; 106 CeedBasis basis_q; 107 108 PetscFunctionBeginUser; 109 // -- Get Pre-requisite things 110 PetscCall(DMGetDimension(projection->dm, &dim)); 111 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q)); 112 113 { // Get newtonian QF context 114 CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op; 115 PetscInt sub_op_index = 0; // will be 0 for the volume op 116 117 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(main_op, &sub_ops)); 118 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 119 } 120 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs)); 121 { // Add the volume integral CeedOperator 122 CeedQFunction qf_rhs_volume; 123 CeedOperator op_rhs_volume; 124 CeedVector q_data; 125 CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 126 CeedBasis basis_diff_flux = NULL; 127 CeedInt q_data_size; 128 129 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 130 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q)); 131 PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q)); 132 PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 133 &elem_restr_qd, &q_data, &q_data_size)); 134 switch (honee->phys->state_var) { 135 case STATEVAR_PRIMITIVE: 136 PetscCallCeed(ceed, 137 CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume)); 138 break; 139 case STATEVAR_CONSERVATIVE: 140 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc, 141 &qf_rhs_volume)); 142 break; 143 case STATEVAR_ENTROPY: 144 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc, 145 &qf_rhs_volume)); 146 break; 147 } 148 149 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx)); 150 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP)); 151 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 152 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 153 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 154 155 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 156 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 157 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 158 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 159 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 160 161 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume)); 162 163 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 164 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 165 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 166 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 167 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 168 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 169 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 170 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 171 } 172 173 { // Add the boundary integral CeedOperator 174 CeedQFunction qf_rhs_boundary; 175 DMLabel face_sets_label; 176 PetscInt num_face_set_values, *face_set_values; 177 CeedInt q_data_size; 178 179 // -- Build RHS operator 180 switch (honee->phys->state_var) { 181 case STATEVAR_PRIMITIVE: 182 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc, 183 &qf_rhs_boundary)); 184 break; 185 case STATEVAR_CONSERVATIVE: 186 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc, 187 &qf_rhs_boundary)); 188 break; 189 case STATEVAR_ENTROPY: 190 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc, 191 &qf_rhs_boundary)); 192 break; 193 } 194 195 PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 196 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx)); 197 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP)); 198 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 199 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 200 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 201 202 PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 203 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 204 for (PetscInt f = 0; f < num_face_set_values; f++) { 205 DMLabel face_orientation_label; 206 PetscInt num_orientations_values, *orientation_values; 207 208 { 209 char *face_orientation_label_name; 210 211 PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 212 PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 213 PetscCall(PetscFree(face_orientation_label_name)); 214 } 215 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 216 for (PetscInt o = 0; o < num_orientations_values; o++) { 217 CeedOperator op_rhs_boundary; 218 CeedBasis basis_q, basis_diff_flux_boundary; 219 CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 220 CeedVector q_data; 221 CeedInt q_data_size; 222 PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 223 224 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 225 PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 226 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 227 &elem_restr_diff_flux_boundary)); 228 PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 229 PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, 230 &q_data_size)); 231 232 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 233 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 234 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 235 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 236 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 237 CEED_VECTOR_ACTIVE)); 238 239 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_boundary)); 240 241 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 242 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 243 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 244 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 245 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 246 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 247 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 248 } 249 PetscCall(PetscFree(orientation_values)); 250 } 251 PetscCall(PetscFree(face_set_values)); 252 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 253 } 254 255 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 256 PetscFunctionReturn(PETSC_SUCCESS); 257 } 258 259 /** 260 @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 261 262 @param[in] honee `Honee` context 263 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 264 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 265 **/ 266 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 267 Ceed ceed = honee->ceed; 268 NodalProjectionData projection = diff_flux_proj->projection; 269 CeedBasis basis_diff_flux, basis_q; 270 CeedElemRestriction elem_restr_diff_flux, elem_restr_qd, elem_restr_q; 271 CeedVector q_data; 272 CeedInt q_data_size; 273 PetscInt dim, num_comp_q; 274 PetscInt height = 0, dm_field = 0; 275 CeedQFunction qf_rhs; 276 CeedQFunctionContext newtonian_qfctx = NULL; 277 278 PetscFunctionBeginUser; 279 PetscCall(DMGetDimension(projection->dm, &dim)); 280 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q)); 281 282 { // Get newtonian QF context 283 CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op; 284 PetscInt sub_op_index = 0; // will be 0 for the volume op 285 286 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(main_op, &sub_ops)); 287 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 288 } 289 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q)); 290 PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q)); 291 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux)); 292 PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux)); 293 PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 294 &elem_restr_qd, &q_data, &q_data_size)); 295 296 switch (honee->phys->state_var) { 297 case STATEVAR_PRIMITIVE: 298 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs)); 299 break; 300 case STATEVAR_CONSERVATIVE: 301 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs)); 302 break; 303 case STATEVAR_ENTROPY: 304 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs)); 305 break; 306 } 307 308 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx)); 309 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP)); 310 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 311 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 312 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 313 314 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 315 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 316 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 317 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 318 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 319 320 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 321 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 322 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 323 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 324 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 325 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 326 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 327 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 328 PetscFunctionReturn(PETSC_SUCCESS); 329 } 330 331 static PetscErrorCode BoundaryIntegralBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 332 HoneeBCStruct honee_bc; 333 334 PetscFunctionBeginUser; 335 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 336 Honee honee = honee_bc->honee; 337 338 switch (honee->phys->state_var) { 339 case STATEVAR_CONSERVATIVE: 340 PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Conserv, BoundaryIntegral_Conserv_loc, honee_bc->qfctx, qf)); 341 break; 342 case STATEVAR_PRIMITIVE: 343 PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Prim, BoundaryIntegral_Prim_loc, honee_bc->qfctx, qf)); 344 break; 345 case STATEVAR_ENTROPY: 346 PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Entropy, BoundaryIntegral_Entropy_loc, honee_bc->qfctx, qf)); 347 break; 348 } 349 PetscFunctionReturn(PETSC_SUCCESS); 350 } 351 352 static PetscErrorCode BoundaryIntegralBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) { 353 HoneeBCStruct honee_bc; 354 355 PetscFunctionBeginUser; 356 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 357 Honee honee = honee_bc->honee; 358 359 switch (honee->phys->state_var) { 360 case STATEVAR_CONSERVATIVE: 361 PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Conserv, BoundaryIntegral_Jacobian_Conserv_loc, honee_bc->qfctx, qf)); 362 break; 363 case STATEVAR_PRIMITIVE: 364 PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Prim, BoundaryIntegral_Jacobian_Prim_loc, honee_bc->qfctx, qf)); 365 break; 366 case STATEVAR_ENTROPY: 367 PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Entropy, BoundaryIntegral_Jacobian_Entropy_loc, honee_bc->qfctx, qf)); 368 break; 369 } 370 PetscFunctionReturn(PETSC_SUCCESS); 371 } 372 373 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx) { 374 SetupContext setup_context; 375 Honee honee = *(Honee *)ctx; 376 CeedInt degree = honee->app_ctx->degree; 377 StabilizationType stab; 378 StateVariable state_var; 379 MPI_Comm comm = honee->comm; 380 Ceed ceed = honee->ceed; 381 PetscBool implicit; 382 PetscBool unit_tests; 383 NewtonianIdealGasContext newtonian_ig_ctx; 384 385 PetscFunctionBeginUser; 386 // Option Defaults 387 const CeedScalar Cv_func[3] = {36, 60, 128}; 388 StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 389 CeedScalar idl_decay_time = -1; 390 PetscCall(PetscNew(&newtonian_ig_ctx)); 391 *newtonian_ig_ctx = (struct NewtonianIdealGasContext_){ 392 .gas = 393 { 394 .cv = 717., 395 .cp = 1004., 396 .lambda = -2. / 3., 397 .mu = 1.8e-5, 398 .k = 0.02638, 399 }, 400 .tau_coeffs = 401 { 402 .Ctau_t = 1.0, 403 .Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1], 404 .Ctau_C = 0.25 / degree, 405 .Ctau_M = 0.25 / degree, 406 .Ctau_E = 0.125, 407 }, 408 .g = {0, 0, 0}, // m/s^2 409 .idl_start = 0, 410 .idl_length = 0, 411 .idl_pressure = reference.pressure, 412 .idl_enable = PETSC_FALSE, 413 }; 414 415 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 416 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 417 (PetscEnum *)&state_var, NULL)); 418 419 // Newtonian fluid properties 420 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, newtonian_ig_ctx->gas.cv, &newtonian_ig_ctx->gas.cv, NULL)); 421 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, newtonian_ig_ctx->gas.cp, &newtonian_ig_ctx->gas.cp, NULL)); 422 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, newtonian_ig_ctx->gas.lambda, 423 &newtonian_ig_ctx->gas.lambda, NULL)); 424 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, newtonian_ig_ctx->gas.mu, &newtonian_ig_ctx->gas.mu, NULL)); 425 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, newtonian_ig_ctx->gas.k, &newtonian_ig_ctx->gas.k, NULL)); 426 427 PetscInt dim = 3; 428 PetscBool given_option = PETSC_FALSE; 429 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 430 PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, newtonian_ig_ctx->g, &dim, &given_option)); 431 if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 432 433 // Stabilization parameters 434 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 435 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_t, 436 &newtonian_ig_ctx->tau_coeffs.Ctau_t, NULL)); 437 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_v, 438 &newtonian_ig_ctx->tau_coeffs.Ctau_v, NULL)); 439 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_C, 440 &newtonian_ig_ctx->tau_coeffs.Ctau_C, NULL)); 441 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_M, 442 &newtonian_ig_ctx->tau_coeffs.Ctau_M, NULL)); 443 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, newtonian_ig_ctx->tau_coeffs.Ctau_E, 444 &newtonian_ig_ctx->tau_coeffs.Ctau_E, NULL)); 445 446 dim = 3; 447 PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 448 PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 449 PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 450 451 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 452 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 453 PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 454 "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 455 456 // IDL Settings 457 { 458 PetscBool idl_enable = (PetscBool)newtonian_ig_ctx->idl_enable; // Need PetscBool variable to read in from PetscOptionsScalar() 459 PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 460 NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 461 PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 462 if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 463 newtonian_ig_ctx->idl_enable = idl_enable; 464 465 PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, newtonian_ig_ctx->idl_start, &newtonian_ig_ctx->idl_start, 466 NULL)); 467 PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, newtonian_ig_ctx->idl_length, 468 &newtonian_ig_ctx->idl_length, NULL)); 469 newtonian_ig_ctx->idl_pressure = reference.pressure; 470 PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, 471 newtonian_ig_ctx->idl_pressure, &newtonian_ig_ctx->idl_pressure, NULL)); 472 } 473 PetscOptionsEnd(); 474 475 // ------------------------------------------------------ 476 // Set up the QFunction context 477 // ------------------------------------------------------ 478 // -- Scale variables to desired units 479 Units units = honee->units; 480 newtonian_ig_ctx->gas.cv *= units->J_per_kg_K; 481 newtonian_ig_ctx->gas.cp *= units->J_per_kg_K; 482 newtonian_ig_ctx->gas.mu *= units->Pascal * units->second; 483 newtonian_ig_ctx->gas.k *= units->W_per_m_K; 484 for (PetscInt i = 0; i < 3; i++) newtonian_ig_ctx->g[i] *= units->m_per_squared_s; 485 reference.pressure *= units->Pascal; 486 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= units->meter / units->second; 487 reference.temperature *= units->Kelvin; 488 489 PetscReal domain_min[3], domain_max[3], domain_size[3]; 490 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 491 for (PetscInt i = 0; i < 3; i++) domain_size[i] = (domain_max[i] - domain_min[i]) * units->meter; 492 493 // -- Solver Settings 494 honee->phys->implicit = implicit; 495 honee->phys->state_var = state_var; 496 497 // -- QFunction Context 498 newtonian_ig_ctx->stabilization = stab; 499 newtonian_ig_ctx->is_implicit = implicit; 500 newtonian_ig_ctx->state_var = state_var; 501 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * units->second); 502 newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 503 504 // -- Setup Context 505 PetscCall(PetscNew(&setup_context)); 506 *setup_context = (struct SetupContext_){ 507 .reference = reference, 508 .newt_ctx = *newtonian_ig_ctx, 509 .lx = domain_size[0], 510 .ly = domain_size[1], 511 .lz = domain_size[2], 512 .time = 0, 513 }; 514 515 CeedQFunctionContext ics_qfctx, newtonian_ig_qfctx; 516 517 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); 518 PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 519 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 520 PetscCallCeed(ceed, 521 CeedQFunctionContextRegisterDouble(ics_qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation")); 522 523 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx)); 524 PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 525 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 526 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 527 "Size of timestep, delta t")); 528 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift", 529 offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 530 "Shift for mass matrix in IJacobian")); 531 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 532 "Current solution time")); 533 534 // Set problem information 535 problem->num_comps_jac_data = 14; 536 if (newtonian_ig_ctx->idl_enable) problem->num_comps_jac_data += 1; 537 problem->compute_exact_solution_error = PETSC_FALSE; 538 problem->print_info = PRINT_NEWTONIAN; 539 problem->num_components = 5; 540 PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 541 static const char *const conserv_component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 542 static const char *const prim_component_names[] = {"Pressure", "VelocityX", "VelocityY", "VelocityZ", "Temperature"}; 543 static const char *const entropy_component_names[] = {"EntropyDensity", "EntropyMomentumX", "EntropyMomentumY", "EntropyMomentumZ", 544 "EntropyTotalEnergy"}; 545 546 switch (state_var) { 547 case STATEVAR_CONSERVATIVE: 548 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Conserv, .qf_loc = ICsNewtonianIG_Conserv_loc}; 549 problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHSFunction_Newtonian, .qf_loc = RHSFunction_Newtonian_loc}; 550 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Conserv, .qf_loc = IFunction_Newtonian_Conserv_loc}; 551 problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Conserv, .qf_loc = IJacobian_Newtonian_Conserv_loc}; 552 for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(conserv_component_names[i], &problem->component_names[i])); 553 break; 554 case STATEVAR_PRIMITIVE: 555 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Prim, .qf_loc = ICsNewtonianIG_Prim_loc}; 556 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Prim, .qf_loc = IFunction_Newtonian_Prim_loc}; 557 problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Prim, .qf_loc = IJacobian_Newtonian_Prim_loc}; 558 for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(prim_component_names[i], &problem->component_names[i])); 559 break; 560 case STATEVAR_ENTROPY: 561 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Entropy, .qf_loc = ICsNewtonianIG_Entropy_loc}; 562 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Entropy, .qf_loc = IFunction_Newtonian_Entropy_loc}; 563 problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Entropy, .qf_loc = IJacobian_Newtonian_Entropy_loc}; 564 for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(entropy_component_names[i], &problem->component_names[i])); 565 break; 566 } 567 // All QFunctions get the same QFunctionContext regardless of state variable 568 problem->ics.qfctx = ics_qfctx; 569 problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; 570 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx)); 571 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx)); 572 573 if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 574 575 PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 4, &honee->diff_flux_proj)); 576 if (honee->diff_flux_proj) { 577 DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 578 NodalProjectionData projection = diff_flux_proj->projection; 579 PetscSection section; 580 581 diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_NS; 582 diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS; 583 PetscCall(DMGetLocalSection(projection->dm, §ion)); 584 switch (honee->diff_flux_proj->method) { 585 case DIV_DIFF_FLUX_PROJ_DIRECT: { 586 PetscCall(PetscSectionSetFieldName(section, 0, "")); 587 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX")); 588 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY")); 589 PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ")); 590 PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy")); 591 } break; 592 case DIV_DIFF_FLUX_PROJ_INDIRECT: { 593 PetscCall(PetscSectionSetFieldName(section, 0, "")); 594 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX")); 595 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY")); 596 PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ")); 597 PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX")); 598 PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY")); 599 PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ")); 600 PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX")); 601 PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY")); 602 PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ")); 603 PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX")); 604 PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY")); 605 PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ")); 606 } break; 607 case DIV_DIFF_FLUX_PROJ_NONE: 608 SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 609 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 610 break; 611 } 612 } 613 614 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 615 BCDefinition bc_def = problem->bc_defs[b]; 616 const char *name; 617 618 PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 619 if (!strcmp(name, "slip")) { 620 PetscCall(SlipBCSetup(bc_def, problem, dm, ctx, newtonian_ig_qfctx)); 621 } else if (!strcmp(name, "freestream")) { 622 PetscCall(FreestreamBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference)); 623 } else if (!strcmp(name, "outflow")) { 624 PetscCall(OutflowBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference)); 625 } else if (!strcmp(name, "inflow")) { 626 HoneeBCStruct honee_bc; 627 628 PetscCall(PetscNew(&honee_bc)); 629 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &honee_bc->qfctx)); 630 honee_bc->honee = honee; 631 honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0; 632 PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 633 634 PetscCall(BCDefinitionSetIFunction(bc_def, BoundaryIntegralBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 635 PetscCall(BCDefinitionSetIJacobian(bc_def, BoundaryIntegralBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp)); 636 } 637 } 638 639 if (unit_tests) PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx->gas)); 640 PetscFunctionReturn(PETSC_SUCCESS); 641 } 642 643 //------------------------------------ 644 // Unit test functions 645 //------------------------------------ 646 647 static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name, 648 PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) { 649 CeedScalar relative_error[5]; // relative error 650 CeedScalar divisor_threshold = 10 * CEED_EPSILON; 651 652 PetscFunctionBeginUser; 653 relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1); 654 relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1); 655 656 CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3])); 657 CeedScalar u_divisor = u_magnitude > divisor_threshold ? u_magnitude : 1; 658 for (int i = 1; i < 4; i++) { 659 relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor; 660 } 661 662 if (fabs(relative_error[0]) >= rtol_0) { 663 printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]); 664 } 665 for (int i = 1; i < 4; i++) { 666 if (fabs(relative_error[i]) >= rtol_u) { 667 printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]); 668 } 669 } 670 if (fabs(relative_error[4]) >= rtol_4) { 671 printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]); 672 } 673 PetscFunctionReturn(PETSC_SUCCESS); 674 } 675 676 // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test 677 static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIGProperties gas, const CeedScalar A0[5], 678 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 679 CeedScalar B0[5], A0_test[5]; 680 char buf[128]; 681 const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 682 683 PetscFunctionBeginUser; 684 const char *A_initial = StateVariables_Initial[state_var_A]; 685 const char *B_initial = StateVariables_Initial[state_var_B]; 686 687 State state_A0 = StateFromQ(gas, A0, state_var_A); 688 StateToQ(gas, state_A0, B0, state_var_B); 689 State state_B0 = StateFromQ(gas, B0, state_var_B); 690 StateToQ(gas, state_B0, A0_test, state_var_A); 691 692 snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial); 693 PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4)); 694 PetscFunctionReturn(PETSC_SUCCESS); 695 } 696 697 // @brief Verify `StateFromQ_fwd` via a finite difference approximation 698 static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIGProperties gas, const CeedScalar A0[5], 699 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 700 CeedScalar eps = 4e-7; // Finite difference step 701 char buf[128]; 702 const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 703 704 PetscFunctionBeginUser; 705 const char *A_initial = StateVariables_Initial[state_var_A]; 706 const char *B_initial = StateVariables_Initial[state_var_B]; 707 State state_0 = StateFromQ(gas, A0, state_var_A); 708 709 for (int i = 0; i < 5; i++) { 710 CeedScalar dB[5] = {0.}, dB_fd[5] = {0.}; 711 { // Calculate dB using State functions 712 CeedScalar dA[5] = {0}; 713 714 dA[i] = A0[i]; 715 State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A); 716 StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B); 717 } 718 719 { // Calculate dB_fd via finite difference approximation 720 CeedScalar A1[5], B0[5], B1[5]; 721 722 for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j]; 723 State state_1 = StateFromQ(gas, A1, state_var_A); 724 StateToQ(gas, state_0, B0, state_var_B); 725 StateToQ(gas, state_1, B1, state_var_B); 726 for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps; 727 } 728 729 snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial); 730 PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4)); 731 } 732 PetscFunctionReturn(PETSC_SUCCESS); 733 } 734 735 // @brief Test the Newtonian State transformation functions, `StateFrom*` 736 static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIGProperties gas) { 737 Units units = honee->units; 738 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin; 739 CeedScalar rtol; 740 741 PetscFunctionBeginUser; 742 const CeedScalar T = 200 * K; 743 const CeedScalar rho = 1.2 * kg / Cube(m); 744 const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas.cv * T; 745 const CeedScalar u_base = 40 * m / sec; 746 const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 747 const CeedScalar e_kinetic = 0.5 * Dot3(u, u); 748 const CeedScalar e_internal = gas.cv * T; 749 const CeedScalar e_total = e_kinetic + e_internal; 750 const CeedScalar gamma = HeatCapacityRatio(gas); 751 const CeedScalar entropy = log(P) - gamma * log(rho); 752 const CeedScalar rho_div_p = rho / P; 753 const CeedScalar Y0[5] = {P, u[0], u[1], u[2], T}; 754 const CeedScalar U0[5] = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total}; 755 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], 756 -rho_div_p}; 757 758 rtol = 20 * CEED_EPSILON; 759 PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 760 PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 761 PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 762 PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol)); 763 PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol)); 764 PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol)); 765 766 rtol = 5e-6; 767 PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 768 PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 769 PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 770 PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol)); 771 PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol)); 772 PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol)); 773 PetscFunctionReturn(PETSC_SUCCESS); 774 } 775