13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 388b783a1SJames Wright // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 588b783a1SJames Wright // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 788b783a1SJames Wright 888b783a1SJames Wright /// @file 988b783a1SJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 1088b783a1SJames Wright 1188b783a1SJames Wright #include "../navierstokes.h" 1288b783a1SJames Wright #include "../qfunctions/setupgeo.h" 1388b783a1SJames Wright #include "../qfunctions/newtonian.h" 1488b783a1SJames Wright 15e334ad8fSJed Brown // Compute relative error |a - b|/|s| 16e334ad8fSJed Brown static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, 17e334ad8fSJed Brown StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 18e334ad8fSJed Brown PetscReal rtol_velocity, PetscReal rtol_temperature) { 196ef2784eSLeila Ghaffari 20e334ad8fSJed Brown PetscFunctionBeginUser; 21e334ad8fSJed Brown StatePrimitive eY; // relative error 22e334ad8fSJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 23e334ad8fSJed Brown PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square( 24e334ad8fSJed Brown sY.velocity[2])); 25e334ad8fSJed Brown for (int j=0; j<3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 26e334ad8fSJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 27a1df05f8SJed Brown if (fabs(eY.pressure) > rtol_pressure) 28a1df05f8SJed Brown printf("%s: pressure error %g\n", name, eY.pressure); 29a1df05f8SJed Brown for (int j=0; j<3; j++) 30a1df05f8SJed Brown if (fabs(eY.velocity[j]) > rtol_velocity) 31e334ad8fSJed Brown printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 32e334ad8fSJed Brown if (fabs(eY.temperature) > rtol_temperature) 33e334ad8fSJed Brown printf("%s: temperature error %g\n", name, eY.temperature); 34e334ad8fSJed Brown PetscFunctionReturn(0); 35e334ad8fSJed Brown } 36e334ad8fSJed Brown 37e334ad8fSJed Brown static PetscErrorCode UnitTests_Newtonian(User user, 38e334ad8fSJed Brown NewtonianIdealGasContext gas) { 396ef2784eSLeila Ghaffari 40e334ad8fSJed Brown Units units = user->units; 41e334ad8fSJed Brown const CeedScalar eps = 1e-6; 42e334ad8fSJed Brown const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, 43e334ad8fSJed Brown Pascal = units->Pascal; 44e334ad8fSJed Brown 45e334ad8fSJed Brown PetscFunctionBeginUser; 46e334ad8fSJed Brown const CeedScalar rho = 1.2 * kg / (m*m*m), u = 40 * m/sec; 47e334ad8fSJed Brown CeedScalar U[5] = {rho, rho*u, rho *u*1.1, rho *u*1.2, 250e3*Pascal + .5*rho *u*u}; 48e334ad8fSJed Brown const CeedScalar x[3] = {.1, .2, .3}; 49e334ad8fSJed Brown State s = StateFromU(gas, U, x); 50e334ad8fSJed Brown for (int i=0; i<8; i++) { 51e334ad8fSJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 52e334ad8fSJed Brown if (i < 5) dU[i] = U[i]; 53e334ad8fSJed Brown else dx[i-5] = x[i-5]; 54e334ad8fSJed Brown State ds = StateFromU_fwd(gas, s, dU, x, dx); 55e334ad8fSJed Brown for (int j=0; j<5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 56e334ad8fSJed Brown for (int j=0; j<3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 57e334ad8fSJed Brown State t = StateFromU(gas, dU, dx); 58e334ad8fSJed Brown StatePrimitive dY; 59e334ad8fSJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 60e334ad8fSJed Brown for (int j=0; j<3; j++) 61e334ad8fSJed Brown dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 62e334ad8fSJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 63e334ad8fSJed Brown char buf[128]; 64e334ad8fSJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 65e334ad8fSJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 66e334ad8fSJed Brown } 67e334ad8fSJed Brown PetscFunctionReturn(0); 68e334ad8fSJed Brown } 69e334ad8fSJed Brown 70a0add3c9SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) { 716ef2784eSLeila Ghaffari 72a0add3c9SJed Brown SetupContext setup_context; 7388b783a1SJames Wright User user = *(User *)ctx; 7488b783a1SJames Wright StabilizationType stab; 7588b783a1SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 7688b783a1SJames Wright PetscBool implicit; 77e334ad8fSJed Brown PetscBool has_curr_time = PETSC_FALSE, unit_tests; 7888b783a1SJames Wright PetscInt ierr; 79841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 80841e4c73SJed Brown CeedQFunctionContext newtonian_ig_context; 8188b783a1SJames Wright 82841e4c73SJed Brown PetscFunctionBeginUser; 83a0add3c9SJed Brown ierr = PetscCalloc1(1, &setup_context); CHKERRQ(ierr); 84841e4c73SJed Brown ierr = PetscCalloc1(1, &newtonian_ig_ctx); CHKERRQ(ierr); 8588b783a1SJames Wright 8688b783a1SJames Wright // ------------------------------------------------------ 8788b783a1SJames Wright // Setup Generic Newtonian IG Problem 8888b783a1SJames Wright // ------------------------------------------------------ 8988b783a1SJames Wright problem->dim = 3; 9088b783a1SJames Wright problem->q_data_size_vol = 10; 91ba6664aeSJames Wright problem->q_data_size_sur = 10; 92*0ec2498eSJames Wright problem->jac_data_size_sur = 11; 9391e5af17SJed Brown problem->setup_vol.qfunction = Setup; 9491e5af17SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 9591e5af17SJed Brown problem->ics.qfunction = ICsNewtonianIG; 9691e5af17SJed Brown problem->ics.qfunction_loc = ICsNewtonianIG_loc; 9791e5af17SJed Brown problem->setup_sur.qfunction = SetupBoundary; 9891e5af17SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 995c677226SJed Brown problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 1005c677226SJed Brown problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 10191e5af17SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Newtonian; 10291e5af17SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_loc; 103e334ad8fSJed Brown problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian; 104e334ad8fSJed Brown problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_loc; 10565dd5cafSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral; 10665dd5cafSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_loc; 10730e9fa81SJames Wright problem->apply_outflow.qfunction = PressureOutflow; 10830e9fa81SJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_loc; 10930e9fa81SJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian; 11030e9fa81SJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_loc; 111a0add3c9SJed Brown problem->bc = NULL; 112a0add3c9SJed Brown problem->bc_ctx = setup_context; 11388b783a1SJames Wright problem->non_zero_time = PETSC_FALSE; 11488b783a1SJames Wright problem->print_info = PRINT_DENSITY_CURRENT; 11588b783a1SJames Wright 11688b783a1SJames Wright // ------------------------------------------------------ 11788b783a1SJames Wright // Create the libCEED context 11888b783a1SJames Wright // ------------------------------------------------------ 11988b783a1SJames Wright CeedScalar cv = 717.; // J/(kg K) 12088b783a1SJames Wright CeedScalar cp = 1004.; // J/(kg K) 12188626eedSJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 12288b783a1SJames Wright CeedScalar lambda = -2./3.; // - 12388626eedSJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 12488b783a1SJames Wright CeedScalar k = 0.02638; // W/(m K) 12588b783a1SJames Wright CeedScalar c_tau = 0.5; // - 12688626eedSJames Wright CeedScalar Ctau_t = 1.0; // - 12788626eedSJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 12888626eedSJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 12988626eedSJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 13088626eedSJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 13188b783a1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 13288b783a1SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 133ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 13488b783a1SJames Wright 13588b783a1SJames Wright // ------------------------------------------------------ 13688b783a1SJames Wright // Create the PETSc context 13788b783a1SJames Wright // ------------------------------------------------------ 13888626eedSJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 13988626eedSJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 14088626eedSJames Wright PetscScalar second = 1; // 1 second in scaled time units 14188b783a1SJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 14288b783a1SJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 14388b783a1SJames Wright 14488b783a1SJames Wright // ------------------------------------------------------ 14588b783a1SJames Wright // Command line Options 14688b783a1SJames Wright // ------------------------------------------------------ 14767490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", 14867490bc6SJeremy L Thompson NULL); 14967490bc6SJeremy L Thompson 15088b783a1SJames Wright // -- Physics 15188b783a1SJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 15288b783a1SJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 15388b783a1SJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 15488b783a1SJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 15588b783a1SJames Wright ierr = PetscOptionsScalar("-lambda", 15688b783a1SJames Wright "Stokes hypothesis second viscosity coefficient", 15788b783a1SJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 15888b783a1SJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 15988b783a1SJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 16088b783a1SJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 16188b783a1SJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 16288b783a1SJames Wright 16388626eedSJames Wright PetscInt dim = problem->dim; 16488626eedSJames Wright ierr = PetscOptionsRealArray("-g", "Gravitational acceleration", 16588626eedSJames Wright NULL, g, &dim, NULL); CHKERRQ(ierr); 16688b783a1SJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 16788b783a1SJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 16888b783a1SJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 16988b783a1SJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 17088b783a1SJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 17188626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant", 17288626eedSJames Wright NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr); 17388626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", 17488626eedSJames Wright NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr); 17588626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", 17688626eedSJames Wright NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr); 17788626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", 17888626eedSJames Wright NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr); 17988626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", 18088626eedSJames Wright NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr); 18188b783a1SJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 18288b783a1SJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 18388b783a1SJames Wright CHKERRQ(ierr); 184e334ad8fSJed Brown ierr = PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", 185e334ad8fSJed Brown NULL, unit_tests=PETSC_FALSE, &unit_tests, NULL); 186e334ad8fSJed Brown CHKERRQ(ierr); 18788b783a1SJames Wright 18888b783a1SJames Wright // -- Units 18988b783a1SJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 19088b783a1SJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 19188b783a1SJames Wright meter = fabs(meter); 19288b783a1SJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 19388b783a1SJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 19488b783a1SJames Wright kilogram = fabs(kilogram); 19588b783a1SJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 19688b783a1SJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 19788b783a1SJames Wright second = fabs(second); 19888b783a1SJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 19988b783a1SJames Wright "1 Kelvin in scaled temperature units", 20088b783a1SJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 20188b783a1SJames Wright Kelvin = fabs(Kelvin); 20288b783a1SJames Wright 20388b783a1SJames Wright // -- Warnings 20488b783a1SJames Wright if (stab == STAB_SUPG && !implicit) { 20588b783a1SJames Wright ierr = PetscPrintf(comm, 20688b783a1SJames Wright "Warning! Use -stab supg only with -implicit\n"); 20788b783a1SJames Wright CHKERRQ(ierr); 20888b783a1SJames Wright } 20967490bc6SJeremy L Thompson PetscOptionsEnd(); 21088b783a1SJames Wright 21188b783a1SJames Wright // ------------------------------------------------------ 21288b783a1SJames Wright // Set up the PETSc context 21388b783a1SJames Wright // ------------------------------------------------------ 21488b783a1SJames Wright // -- Define derived units 21588b783a1SJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 21688b783a1SJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 21788b783a1SJames Wright m_per_squared_s = meter / PetscSqr(second); 21888b783a1SJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 21988b783a1SJames Wright 22088b783a1SJames Wright user->units->meter = meter; 22188b783a1SJames Wright user->units->kilogram = kilogram; 22288b783a1SJames Wright user->units->second = second; 22388b783a1SJames Wright user->units->Kelvin = Kelvin; 22488b783a1SJames Wright user->units->Pascal = Pascal; 22588b783a1SJames Wright user->units->J_per_kg_K = J_per_kg_K; 22688b783a1SJames Wright user->units->m_per_squared_s = m_per_squared_s; 22788b783a1SJames Wright user->units->W_per_m_K = W_per_m_K; 22888b783a1SJames Wright 22988b783a1SJames Wright // ------------------------------------------------------ 23088b783a1SJames Wright // Set up the libCEED context 23188b783a1SJames Wright // ------------------------------------------------------ 23288b783a1SJames Wright // -- Scale variables to desired units 23388b783a1SJames Wright cv *= J_per_kg_K; 23488b783a1SJames Wright cp *= J_per_kg_K; 23588b783a1SJames Wright mu *= Pascal * second; 23688b783a1SJames Wright k *= W_per_m_K; 237ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] *= meter; 238ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) g[i] *= m_per_squared_s; 23988b783a1SJames Wright problem->dm_scale = meter; 24088b783a1SJames Wright 24188b783a1SJames Wright // -- Setup Context 24288b783a1SJames Wright setup_context->cv = cv; 24388b783a1SJames Wright setup_context->cp = cp; 24488b783a1SJames Wright setup_context->lx = domain_size[0]; 24588b783a1SJames Wright setup_context->ly = domain_size[1]; 24688b783a1SJames Wright setup_context->lz = domain_size[2]; 24788b783a1SJames Wright setup_context->time = 0; 24888626eedSJames Wright ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr); 24988b783a1SJames Wright 25088b783a1SJames Wright // -- Solver Settings 25188b783a1SJames Wright user->phys->stab = stab; 25288b783a1SJames Wright user->phys->implicit = implicit; 25388b783a1SJames Wright user->phys->has_curr_time = has_curr_time; 25488b783a1SJames Wright 25588b783a1SJames Wright // -- QFunction Context 256841e4c73SJed Brown newtonian_ig_ctx->lambda = lambda; 257841e4c73SJed Brown newtonian_ig_ctx->mu = mu; 258841e4c73SJed Brown newtonian_ig_ctx->k = k; 259841e4c73SJed Brown newtonian_ig_ctx->cv = cv; 260841e4c73SJed Brown newtonian_ig_ctx->cp = cp; 261841e4c73SJed Brown newtonian_ig_ctx->c_tau = c_tau; 262841e4c73SJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 263841e4c73SJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 264841e4c73SJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 265841e4c73SJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 266841e4c73SJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 267841e4c73SJed Brown newtonian_ig_ctx->stabilization = stab; 26865dd5cafSJames Wright newtonian_ig_ctx->is_implicit = implicit; 269841e4c73SJed Brown ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 27088b783a1SJames Wright 271841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 272841e4c73SJed Brown CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, 273841e4c73SJed Brown CEED_USE_POINTER, sizeof(*setup_context), setup_context); 27417ce10faSJeremy L Thompson CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, 27517ce10faSJeremy L Thompson CEED_MEM_HOST, 27617ce10faSJeremy L Thompson FreeContextPetsc); 277841e4c73SJed Brown CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, 278841e4c73SJed Brown "evaluation time", 279841e4c73SJed Brown (char *)&setup_context->time - (char *)setup_context, 1, "Time of evaluation"); 280841e4c73SJed Brown 281841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 282841e4c73SJed Brown CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, 283841e4c73SJed Brown CEED_USE_POINTER, 284841e4c73SJed Brown sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 285841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, 286841e4c73SJed Brown FreeContextPetsc); 287841e4c73SJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", 288841e4c73SJed Brown offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t"); 289e334ad8fSJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 290e334ad8fSJed Brown offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 291e334ad8fSJed Brown "Shift for mass matrix in IJacobian"); 292841e4c73SJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 293841e4c73SJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 294841e4c73SJed Brown &problem->apply_vol_ifunction.qfunction_context); 295e334ad8fSJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 296e334ad8fSJed Brown &problem->apply_vol_ijacobian.qfunction_context); 29765dd5cafSJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 29865dd5cafSJames Wright &problem->apply_inflow.qfunction_context); 29930e9fa81SJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 30030e9fa81SJames Wright &problem->apply_outflow.qfunction_context); 30130e9fa81SJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 30230e9fa81SJames Wright &problem->apply_outflow_jacobian.qfunction_context); 303e334ad8fSJed Brown 304e334ad8fSJed Brown if (unit_tests) { 305e334ad8fSJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 306e334ad8fSJed Brown } 30788b783a1SJames Wright PetscFunctionReturn(0); 30888b783a1SJames Wright } 30988b783a1SJames Wright 310841e4c73SJed Brown PetscErrorCode PRINT_DENSITY_CURRENT(ProblemData *problem, 311841e4c73SJed Brown AppCtx app_ctx) { 312841e4c73SJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 313841e4c73SJed Brown PetscErrorCode ierr; 314841e4c73SJed Brown NewtonianIdealGasContext newtonian_ctx; 315841e4c73SJed Brown 31688b783a1SJames Wright PetscFunctionBeginUser; 317841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 318841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ctx); 319841e4c73SJed Brown ierr = PetscPrintf(comm, 320841e4c73SJed Brown " Problem:\n" 321841e4c73SJed Brown " Problem Name : %s\n" 322841e4c73SJed Brown " Stabilization : %s\n", 323841e4c73SJed Brown app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]); 324841e4c73SJed Brown CHKERRQ(ierr); 325841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 326841e4c73SJed Brown &newtonian_ctx); 32788b783a1SJames Wright PetscFunctionReturn(0); 32888b783a1SJames Wright } 329