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(User user, NewtonianIdealGasContext gas) { 107 Units units = user->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(User user, CeedOperator *op_mass) { 155 Ceed ceed = user->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 qf_ctx = 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(user->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], &qf_ctx)); 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, qf_ctx)); 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, user->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(&qf_ctx)); 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] user `User` object 211 @param[in] ceed_data `CeedData` object 212 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 213 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 214 **/ 215 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj, 216 CeedOperator *op_rhs) { 217 Ceed ceed = user->ceed; 218 NodalProjectionData projection = diff_flux_proj->projection; 219 CeedBasis basis_diff_flux; 220 CeedElemRestriction elem_restr_diff_flux_volume; 221 CeedInt num_comp_q; 222 PetscInt dim, label_value = 0; 223 DMLabel domain_label = NULL; 224 CeedQFunctionContext newtonian_qf_ctx = NULL; 225 226 PetscFunctionBeginUser; 227 // -- Get Pre-requisite things 228 PetscCall(DMGetDimension(projection->dm, &dim)); 229 PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q)); 230 231 { // Get elem_restr_diff_flux and basis_diff_flux 232 CeedOperator *sub_ops; 233 CeedOperatorField op_field; 234 PetscInt sub_op_index = 0; // will be 0 for the volume op 235 236 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops)); 237 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "div F_diff", &op_field)); 238 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_diff_flux_volume)); 239 PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_diff_flux)); 240 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qf_ctx)); 241 } 242 PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 243 { // Add the volume integral CeedOperator 244 CeedQFunction qf_rhs_volume; 245 CeedOperator op_rhs_volume; 246 CeedVector q_data; 247 CeedElemRestriction elem_restr_qd; 248 CeedInt q_data_size; 249 250 PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, 251 &elem_restr_qd, &q_data, &q_data_size)); 252 switch (user->phys->state_var) { 253 case STATEVAR_PRIMITIVE: 254 PetscCallCeed(ceed, 255 CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume)); 256 break; 257 case STATEVAR_CONSERVATIVE: 258 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc, 259 &qf_rhs_volume)); 260 break; 261 case STATEVAR_ENTROPY: 262 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc, 263 &qf_rhs_volume)); 264 break; 265 } 266 267 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qf_ctx)); 268 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP)); 269 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 270 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 271 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 272 273 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 274 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 275 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 276 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 277 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 278 279 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 280 281 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 282 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 283 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 284 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 285 } 286 287 { // Add the boundary integral CeedOperator 288 CeedQFunction qf_rhs_boundary; 289 DMLabel face_sets_label; 290 PetscInt num_face_set_values, *face_set_values; 291 CeedInt q_data_size; 292 293 // -- Build RHS operator 294 switch (user->phys->state_var) { 295 case STATEVAR_PRIMITIVE: 296 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc, 297 &qf_rhs_boundary)); 298 break; 299 case STATEVAR_CONSERVATIVE: 300 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc, 301 &qf_rhs_boundary)); 302 break; 303 case STATEVAR_ENTROPY: 304 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc, 305 &qf_rhs_boundary)); 306 break; 307 } 308 309 PetscCall(QDataBoundaryGradientGetNumComponents(user->dm, &q_data_size)); 310 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qf_ctx)); 311 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP)); 312 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 313 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 314 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 315 316 PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 317 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 318 for (PetscInt f = 0; f < num_face_set_values; f++) { 319 DMLabel face_orientation_label; 320 PetscInt num_orientations_values, *orientation_values; 321 322 { 323 char *face_orientation_label_name; 324 325 PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 326 PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 327 PetscCall(PetscFree(face_orientation_label_name)); 328 } 329 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 330 for (PetscInt o = 0; o < num_orientations_values; o++) { 331 CeedOperator op_rhs_boundary; 332 CeedBasis basis_q, basis_diff_flux_boundary; 333 CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 334 CeedVector q_data; 335 CeedInt q_data_size; 336 PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 337 338 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, user->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 339 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 340 &elem_restr_diff_flux_boundary)); 341 PetscCall(QDataBoundaryGradientGet(ceed, user->dm, face_orientation_label, orientation, ceed_data->x_coord, &elem_restr_qdata, &q_data, 342 &q_data_size)); 343 PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, user->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 344 PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 345 346 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 347 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 348 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 349 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 350 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 351 CEED_VECTOR_ACTIVE)); 352 353 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 354 355 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 356 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 357 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 358 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 359 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 360 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 361 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 362 } 363 PetscCall(PetscFree(orientation_values)); 364 } 365 PetscCall(PetscFree(face_set_values)); 366 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 367 } 368 369 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qf_ctx)); 370 PetscFunctionReturn(PETSC_SUCCESS); 371 } 372 373 /** 374 @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 375 376 @param[in] user `User` object 377 @param[in] ceed_data `CeedData` object 378 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 379 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 380 **/ 381 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj, 382 CeedOperator *op_rhs) { 383 Ceed ceed = user->ceed; 384 NodalProjectionData projection = diff_flux_proj->projection; 385 CeedBasis basis_diff_flux; 386 CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 387 CeedVector q_data; 388 CeedInt num_comp_q, q_data_size; 389 PetscInt dim; 390 PetscInt label_value = 0, height = 0, dm_field = 0; 391 DMLabel domain_label = NULL; 392 CeedQFunction qf_rhs; 393 CeedQFunctionContext newtonian_qf_ctx = NULL; 394 395 PetscFunctionBeginUser; 396 PetscCall(DMGetDimension(projection->dm, &dim)); 397 PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q)); 398 399 { // Get elem_restr_div_diff_flux 400 CeedOperator *sub_ops; 401 PetscInt sub_op_index = 0; // will be 0 for the volume op 402 403 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops)); 404 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qf_ctx)); 405 } 406 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 407 PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 408 PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &elem_restr_qd, 409 &q_data, &q_data_size)); 410 411 switch (user->phys->state_var) { 412 case STATEVAR_PRIMITIVE: 413 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs)); 414 break; 415 case STATEVAR_CONSERVATIVE: 416 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs)); 417 break; 418 case STATEVAR_ENTROPY: 419 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs)); 420 break; 421 } 422 423 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qf_ctx)); 424 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP)); 425 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 426 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 427 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 428 429 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 430 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 431 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 432 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 433 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 434 435 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 436 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qf_ctx)); 437 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 438 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 439 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 440 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 441 PetscFunctionReturn(PETSC_SUCCESS); 442 } 443 444 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 445 SetupContext setup_context; 446 User user = *(User *)ctx; 447 CeedInt degree = user->app_ctx->degree; 448 StabilizationType stab; 449 StateVariable state_var; 450 MPI_Comm comm = user->comm; 451 Ceed ceed = user->ceed; 452 PetscBool implicit; 453 PetscBool unit_tests; 454 NewtonianIdealGasContext newtonian_ig_ctx; 455 CeedQFunctionContext newtonian_ig_qfctx; 456 457 PetscFunctionBeginUser; 458 PetscCall(PetscCalloc1(1, &setup_context)); 459 PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 460 461 // ------------------------------------------------------ 462 // Setup Generic Newtonian IG Problem 463 // ------------------------------------------------------ 464 problem->jac_data_size_vol = 14; 465 problem->jac_data_size_sur = 11; 466 problem->compute_exact_solution_error = PETSC_FALSE; 467 problem->print_info = PRINT_NEWTONIAN; 468 469 PetscCall(DivDiffFluxProjectionCreate(user, 4, &user->diff_flux_proj)); 470 if (user->diff_flux_proj) { 471 DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj; 472 NodalProjectionData projection = diff_flux_proj->projection; 473 474 diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_NS; 475 diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS; 476 477 switch (user->diff_flux_proj->method) { 478 case DIV_DIFF_FLUX_PROJ_DIRECT: { 479 PetscSection section; 480 481 PetscCall(DMGetLocalSection(projection->dm, §ion)); 482 PetscCall(PetscSectionSetFieldName(section, 0, "")); 483 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX")); 484 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY")); 485 PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ")); 486 PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy")); 487 } break; 488 case DIV_DIFF_FLUX_PROJ_INDIRECT: { 489 PetscSection section; 490 491 PetscCall(DMGetLocalSection(projection->dm, §ion)); 492 PetscCall(PetscSectionSetFieldName(section, 0, "")); 493 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX")); 494 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY")); 495 PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ")); 496 PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX")); 497 PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY")); 498 PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ")); 499 PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX")); 500 PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY")); 501 PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ")); 502 PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX")); 503 PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY")); 504 PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ")); 505 } break; 506 case DIV_DIFF_FLUX_PROJ_NONE: 507 SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 508 DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]); 509 break; 510 } 511 } 512 513 // ------------------------------------------------------ 514 // Create the libCEED context 515 // ------------------------------------------------------ 516 CeedScalar cv = 717.; // J/(kg K) 517 CeedScalar cp = 1004.; // J/(kg K) 518 CeedScalar g[3] = {0, 0, 0}; // m/s^2 519 CeedScalar lambda = -2. / 3.; // - 520 CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 521 CeedScalar k = 0.02638; // W/(m K) 522 CeedScalar c_tau = 0.5 / degree; // - 523 CeedScalar Ctau_t = 1.0; // - 524 CeedScalar Cv_func[3] = {36, 60, 128}; 525 CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 526 CeedScalar Ctau_C = 0.25 / degree; 527 CeedScalar Ctau_M = 0.25 / degree; 528 CeedScalar Ctau_E = 0.125; 529 PetscReal domain_min[3], domain_max[3], domain_size[3]; 530 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 531 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 532 533 StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 534 CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure; 535 PetscBool idl_enable = PETSC_FALSE; 536 537 // ------------------------------------------------------ 538 // Create the PETSc context 539 // ------------------------------------------------------ 540 PetscScalar meter = 1; // 1 meter in scaled length units 541 PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 542 PetscScalar second = 1; // 1 second in scaled time units 543 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 544 PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 545 546 // ------------------------------------------------------ 547 // Command line Options 548 // ------------------------------------------------------ 549 PetscBool given_option = PETSC_FALSE; 550 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 551 // -- Conservative vs Primitive variables 552 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 553 (PetscEnum *)&state_var, NULL)); 554 555 switch (state_var) { 556 case STATEVAR_CONSERVATIVE: 557 problem->ics.qf_func_ptr = ICsNewtonianIG_Conserv; 558 problem->ics.qf_loc = ICsNewtonianIG_Conserv_loc; 559 problem->apply_vol_rhs.qf_func_ptr = RHSFunction_Newtonian; 560 problem->apply_vol_rhs.qf_loc = RHSFunction_Newtonian_loc; 561 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Conserv; 562 problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Conserv_loc; 563 problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Conserv; 564 problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Conserv_loc; 565 problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Conserv; 566 problem->apply_inflow.qf_loc = BoundaryIntegral_Conserv_loc; 567 problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Conserv; 568 problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Conserv_loc; 569 break; 570 case STATEVAR_PRIMITIVE: 571 problem->ics.qf_func_ptr = ICsNewtonianIG_Prim; 572 problem->ics.qf_loc = ICsNewtonianIG_Prim_loc; 573 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Prim; 574 problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Prim_loc; 575 problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Prim; 576 problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Prim_loc; 577 problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Prim; 578 problem->apply_inflow.qf_loc = BoundaryIntegral_Prim_loc; 579 problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Prim; 580 problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Prim_loc; 581 break; 582 case STATEVAR_ENTROPY: 583 problem->ics.qf_func_ptr = ICsNewtonianIG_Entropy; 584 problem->ics.qf_loc = ICsNewtonianIG_Entropy_loc; 585 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Entropy; 586 problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Entropy_loc; 587 problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Entropy; 588 problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Entropy_loc; 589 problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Entropy; 590 problem->apply_inflow.qf_loc = BoundaryIntegral_Entropy_loc; 591 problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Entropy; 592 problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Entropy_loc; 593 break; 594 } 595 596 // -- Physics 597 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 598 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 599 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 600 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 601 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 602 603 PetscInt dim = 3; 604 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 605 PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 606 if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 607 608 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 609 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 610 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 611 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 612 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 613 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 614 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 615 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 616 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 617 618 dim = 3; 619 PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 620 PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 621 PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 622 623 // -- Units 624 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 625 meter = fabs(meter); 626 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 627 kilogram = fabs(kilogram); 628 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 629 second = fabs(second); 630 PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 631 Kelvin = fabs(Kelvin); 632 633 // -- Warnings 634 PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 635 "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 636 PetscCheck(!(user->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE && !implicit), comm, PETSC_ERR_SUP, 637 "Projection of divergence of diffusive flux is not implemented for explicit timestepping"); 638 639 PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 640 NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 641 PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 642 if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 643 if (idl_enable) problem->jac_data_size_vol++; 644 PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 645 PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 646 idl_pressure = reference.pressure; 647 PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure, 648 &idl_pressure, NULL)); 649 PetscOptionsEnd(); 650 651 if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 652 653 // ------------------------------------------------------ 654 // Set up the PETSc context 655 // ------------------------------------------------------ 656 // -- Define derived units 657 Pascal = kilogram / (meter * PetscSqr(second)); 658 J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 659 m_per_squared_s = meter / PetscSqr(second); 660 W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 661 662 user->units->meter = meter; 663 user->units->kilogram = kilogram; 664 user->units->second = second; 665 user->units->Kelvin = Kelvin; 666 user->units->Pascal = Pascal; 667 user->units->J_per_kg_K = J_per_kg_K; 668 user->units->m_per_squared_s = m_per_squared_s; 669 user->units->W_per_m_K = W_per_m_K; 670 671 // ------------------------------------------------------ 672 // Set up the libCEED context 673 // ------------------------------------------------------ 674 // -- Scale variables to desired units 675 cv *= J_per_kg_K; 676 cp *= J_per_kg_K; 677 mu *= Pascal * second; 678 k *= W_per_m_K; 679 for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 680 for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 681 reference.pressure *= Pascal; 682 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 683 reference.temperature *= Kelvin; 684 685 // -- Solver Settings 686 user->phys->implicit = implicit; 687 user->phys->state_var = state_var; 688 689 // -- QFunction Context 690 newtonian_ig_ctx->lambda = lambda; 691 newtonian_ig_ctx->mu = mu; 692 newtonian_ig_ctx->k = k; 693 newtonian_ig_ctx->cv = cv; 694 newtonian_ig_ctx->cp = cp; 695 newtonian_ig_ctx->c_tau = c_tau; 696 newtonian_ig_ctx->Ctau_t = Ctau_t; 697 newtonian_ig_ctx->Ctau_v = Ctau_v; 698 newtonian_ig_ctx->Ctau_C = Ctau_C; 699 newtonian_ig_ctx->Ctau_M = Ctau_M; 700 newtonian_ig_ctx->Ctau_E = Ctau_E; 701 newtonian_ig_ctx->stabilization = stab; 702 newtonian_ig_ctx->is_implicit = implicit; 703 newtonian_ig_ctx->state_var = state_var; 704 newtonian_ig_ctx->idl_enable = idl_enable; 705 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 706 newtonian_ig_ctx->idl_start = idl_start * meter; 707 newtonian_ig_ctx->idl_length = idl_length * meter; 708 newtonian_ig_ctx->idl_pressure = idl_pressure; 709 newtonian_ig_ctx->divFdiff_method = user->app_ctx->divFdiffproj_method; 710 PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 711 712 // -- Setup Context 713 setup_context->reference = reference; 714 setup_context->gas = *newtonian_ig_ctx; 715 setup_context->lx = domain_size[0]; 716 setup_context->ly = domain_size[1]; 717 setup_context->lz = domain_size[2]; 718 setup_context->time = 0; 719 720 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfctx)); 721 PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 722 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 723 PetscCallCeed( 724 ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation")); 725 726 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_qfctx)); 727 PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 728 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 729 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 730 "Size of timestep, delta t")); 731 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift", 732 offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 733 "Shift for mass matrix in IJacobian")); 734 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 735 "Current solution time")); 736 737 problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; 738 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx)); 739 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx)); 740 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow.qfctx)); 741 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow_jacobian.qfctx)); 742 743 if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 744 if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 745 if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_qfctx)); 746 747 if (unit_tests) { 748 PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 749 } 750 PetscFunctionReturn(PETSC_SUCCESS); 751 } 752 753 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) { 754 MPI_Comm comm = user->comm; 755 Ceed ceed = user->ceed; 756 NewtonianIdealGasContext newtonian_ctx; 757 758 PetscFunctionBeginUser; 759 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx)); 760 PetscCall(PetscPrintf(comm, 761 " Problem:\n" 762 " Problem Name : %s\n" 763 " Stabilization : %s\n", 764 app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 765 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx)); 766 PetscFunctionReturn(PETSC_SUCCESS); 767 } 768