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