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