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