1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utility functions for setting up problems using the Newtonian Qfunction 10 11 #include "../qfunctions/newtonian.h" 12 13 #include <ceed.h> 14 #include <petscdm.h> 15 16 #include "../navierstokes.h" 17 #include "../qfunctions/setupgeo.h" 18 19 // For use with PetscOptionsEnum 20 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 21 22 // Compute relative error |a - b|/|s| 23 static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 24 PetscReal rtol_velocity, PetscReal rtol_temperature) { 25 StatePrimitive eY; // relative error 26 27 PetscFunctionBeginUser; 28 eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 29 PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 30 for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 31 eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 32 if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 33 for (int j = 0; j < 3; j++) { 34 if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 35 } 36 if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 37 PetscFunctionReturn(PETSC_SUCCESS); 38 } 39 40 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 41 Units units = user->units; 42 const CeedScalar eps = 1e-6; 43 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 44 PetscFunctionBeginUser; 45 const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 46 CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 47 State s = StateFromU(gas, U); 48 for (int i = 0; i < 8; i++) { 49 CeedScalar dU[5] = {0}; 50 if (i < 5) dU[i] = U[i]; 51 State ds = StateFromU_fwd(gas, s, dU); 52 for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 53 State t = StateFromU(gas, dU); 54 StatePrimitive dY; 55 dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 56 for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 57 dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 58 char buf[128]; 59 snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 60 PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 61 } 62 PetscFunctionReturn(PETSC_SUCCESS); 63 } 64 65 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 66 // 67 // Only used for SUPG stabilization 68 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) { 69 Ceed ceed = user->ceed; 70 CeedInt num_comp_q, q_data_size; 71 CeedQFunction qf_mass; 72 CeedElemRestriction elem_restr_q, elem_restr_qd_i; 73 CeedBasis basis_q; 74 CeedVector q_data; 75 CeedQFunctionContext qf_ctx = NULL; 76 PetscInt dim = 3; 77 78 PetscFunctionBeginUser; 79 { // Get restriction and basis from the RHS function 80 CeedOperator *sub_ops; 81 CeedOperatorField field; 82 PetscInt sub_op_index = 0; // will be 0 for the volume op 83 84 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 85 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 86 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 87 PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 88 89 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 90 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 91 PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 92 93 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 94 } 95 96 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 97 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 98 99 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 100 101 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 102 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 103 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 104 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 105 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 106 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 107 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 108 109 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 110 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 111 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 112 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 113 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 114 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 115 116 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 117 PetscFunctionReturn(PETSC_SUCCESS); 118 } 119 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 120 SetupContext setup_context; 121 User user = *(User *)ctx; 122 CeedInt degree = user->app_ctx->degree; 123 StabilizationType stab; 124 StateVariable state_var; 125 MPI_Comm comm = user->comm; 126 Ceed ceed = user->ceed; 127 PetscBool implicit; 128 PetscBool unit_tests; 129 NewtonianIdealGasContext newtonian_ig_ctx; 130 CeedQFunctionContext newtonian_ig_context; 131 132 PetscFunctionBeginUser; 133 PetscCall(PetscCalloc1(1, &setup_context)); 134 PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 135 136 // ------------------------------------------------------ 137 // Setup Generic Newtonian IG Problem 138 // ------------------------------------------------------ 139 problem->dim = 3; 140 problem->q_data_size_vol = 10; 141 problem->q_data_size_sur = 10; 142 problem->jac_data_size_sur = 11; 143 problem->setup_vol.qfunction = Setup; 144 problem->setup_vol.qfunction_loc = Setup_loc; 145 problem->setup_sur.qfunction = SetupBoundary; 146 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 147 problem->non_zero_time = PETSC_FALSE; 148 problem->print_info = PRINT_NEWTONIAN; 149 problem->uses_newtonian = PETSC_TRUE; 150 151 // ------------------------------------------------------ 152 // Create the libCEED context 153 // ------------------------------------------------------ 154 CeedScalar cv = 717.; // J/(kg K) 155 CeedScalar cp = 1004.; // J/(kg K) 156 CeedScalar g[3] = {0, 0, 0}; // m/s^2 157 CeedScalar lambda = -2. / 3.; // - 158 CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 159 CeedScalar k = 0.02638; // W/(m K) 160 CeedScalar c_tau = 0.5 / degree; // - 161 CeedScalar Ctau_t = 1.0; // - 162 CeedScalar Cv_func[3] = {36, 60, 128}; 163 CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 164 CeedScalar Ctau_C = 0.25 / degree; 165 CeedScalar Ctau_M = 0.25 / degree; 166 CeedScalar Ctau_E = 0.125; 167 PetscReal domain_min[3], domain_max[3], domain_size[3]; 168 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 169 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 170 171 StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 172 CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0; 173 PetscBool idl_enable = PETSC_FALSE; 174 175 // ------------------------------------------------------ 176 // Create the PETSc context 177 // ------------------------------------------------------ 178 PetscScalar meter = 1; // 1 meter in scaled length units 179 PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 180 PetscScalar second = 1; // 1 second in scaled time units 181 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 182 PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 183 184 // ------------------------------------------------------ 185 // Command line Options 186 // ------------------------------------------------------ 187 PetscBool given_option = PETSC_FALSE; 188 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 189 // -- Conservative vs Primitive variables 190 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 191 (PetscEnum *)&state_var, NULL)); 192 193 switch (state_var) { 194 case STATEVAR_CONSERVATIVE: 195 problem->ics.qfunction = ICsNewtonianIG_Conserv; 196 problem->ics.qfunction_loc = ICsNewtonianIG_Conserv_loc; 197 problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 198 problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 199 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 200 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 201 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 202 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 203 problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 204 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 205 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 206 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 207 break; 208 209 case STATEVAR_PRIMITIVE: 210 problem->ics.qfunction = ICsNewtonianIG_Prim; 211 problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 212 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 213 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 214 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 215 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 216 problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 217 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 218 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 219 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 220 break; 221 } 222 223 // -- Physics 224 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 225 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 226 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 227 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 228 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 229 230 PetscInt dim = problem->dim; 231 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 232 PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 233 dim = problem->dim; 234 PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 235 if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 236 237 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 238 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 239 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 240 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 241 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 242 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 243 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 244 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 245 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 246 247 dim = 3; 248 PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 249 PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 250 PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 251 252 // -- Units 253 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 254 meter = fabs(meter); 255 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 256 kilogram = fabs(kilogram); 257 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 258 second = fabs(second); 259 PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 260 Kelvin = fabs(Kelvin); 261 262 // -- Warnings 263 PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 264 "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 265 266 PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 267 NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 268 PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 269 if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 270 PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 271 PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 272 PetscOptionsEnd(); 273 274 if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 275 276 // ------------------------------------------------------ 277 // Set up the PETSc context 278 // ------------------------------------------------------ 279 // -- Define derived units 280 Pascal = kilogram / (meter * PetscSqr(second)); 281 J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 282 m_per_squared_s = meter / PetscSqr(second); 283 W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 284 285 user->units->meter = meter; 286 user->units->kilogram = kilogram; 287 user->units->second = second; 288 user->units->Kelvin = Kelvin; 289 user->units->Pascal = Pascal; 290 user->units->J_per_kg_K = J_per_kg_K; 291 user->units->m_per_squared_s = m_per_squared_s; 292 user->units->W_per_m_K = W_per_m_K; 293 294 // ------------------------------------------------------ 295 // Set up the libCEED context 296 // ------------------------------------------------------ 297 // -- Scale variables to desired units 298 cv *= J_per_kg_K; 299 cp *= J_per_kg_K; 300 mu *= Pascal * second; 301 k *= W_per_m_K; 302 for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 303 for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 304 reference.pressure *= Pascal; 305 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 306 reference.temperature *= Kelvin; 307 problem->dm_scale = meter; 308 309 // -- Solver Settings 310 user->phys->implicit = implicit; 311 user->phys->state_var = state_var; 312 313 // -- QFunction Context 314 newtonian_ig_ctx->lambda = lambda; 315 newtonian_ig_ctx->mu = mu; 316 newtonian_ig_ctx->k = k; 317 newtonian_ig_ctx->cv = cv; 318 newtonian_ig_ctx->cp = cp; 319 newtonian_ig_ctx->c_tau = c_tau; 320 newtonian_ig_ctx->Ctau_t = Ctau_t; 321 newtonian_ig_ctx->Ctau_v = Ctau_v; 322 newtonian_ig_ctx->Ctau_C = Ctau_C; 323 newtonian_ig_ctx->Ctau_M = Ctau_M; 324 newtonian_ig_ctx->Ctau_E = Ctau_E; 325 newtonian_ig_ctx->P0 = reference.pressure; 326 newtonian_ig_ctx->stabilization = stab; 327 newtonian_ig_ctx->P0 = reference.pressure; 328 newtonian_ig_ctx->is_implicit = implicit; 329 newtonian_ig_ctx->state_var = state_var; 330 newtonian_ig_ctx->idl_enable = idl_enable; 331 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 332 newtonian_ig_ctx->idl_start = idl_start * meter; 333 newtonian_ig_ctx->idl_length = idl_length * meter; 334 PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 335 336 // -- Setup Context 337 setup_context->reference = reference; 338 setup_context->gas = *newtonian_ig_ctx; 339 setup_context->lx = domain_size[0]; 340 setup_context->ly = domain_size[1]; 341 setup_context->lz = domain_size[2]; 342 setup_context->time = 0; 343 344 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 345 PetscCallCeed(ceed, 346 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 347 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 348 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1, 349 "Time of evaluation")); 350 351 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context)); 352 PetscCallCeed(ceed, 353 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 354 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc)); 355 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 356 "Size of timestep, delta t")); 357 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 358 offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 359 "Shift for mass matrix in IJacobian")); 360 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 361 "Current solution time")); 362 363 problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 364 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context)); 365 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context)); 366 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context)); 367 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context)); 368 369 if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 370 if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 371 if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context)); 372 373 if (unit_tests) { 374 PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 375 } 376 PetscFunctionReturn(PETSC_SUCCESS); 377 } 378 379 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) { 380 MPI_Comm comm = user->comm; 381 Ceed ceed = user->ceed; 382 NewtonianIdealGasContext newtonian_ctx; 383 384 PetscFunctionBeginUser; 385 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx)); 386 PetscCall(PetscPrintf(comm, 387 " Problem:\n" 388 " Problem Name : %s\n" 389 " Stabilization : %s\n", 390 app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 391 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx)); 392 PetscFunctionReturn(PETSC_SUCCESS); 393 } 394