1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 33a8779fbSJames Wright // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 53a8779fbSJames Wright // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73a8779fbSJames Wright 83a8779fbSJames Wright /// @file 93a8779fbSJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 103a8779fbSJames Wright 113a8779fbSJames Wright #include "../navierstokes.h" 123a8779fbSJames Wright #include "../qfunctions/setupgeo.h" 133a8779fbSJames Wright #include "../qfunctions/newtonian.h" 143a8779fbSJames Wright 15f0b65372SJed Brown // Compute relative error |a - b|/|s| 16f0b65372SJed Brown static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, 17f0b65372SJed Brown StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 18f0b65372SJed Brown PetscReal rtol_velocity, PetscReal rtol_temperature) { 196dd99bedSLeila Ghaffari 20f0b65372SJed Brown PetscFunctionBeginUser; 21f0b65372SJed Brown StatePrimitive eY; // relative error 22f0b65372SJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 23f0b65372SJed Brown PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square( 24f0b65372SJed Brown sY.velocity[2])); 25f0b65372SJed Brown for (int j=0; j<3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 26f0b65372SJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 27edd152dcSJed Brown if (fabs(eY.pressure) > rtol_pressure) 28edd152dcSJed Brown printf("%s: pressure error %g\n", name, eY.pressure); 29edd152dcSJed Brown for (int j=0; j<3; j++) 30edd152dcSJed Brown if (fabs(eY.velocity[j]) > rtol_velocity) 31f0b65372SJed Brown printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 32f0b65372SJed Brown if (fabs(eY.temperature) > rtol_temperature) 33f0b65372SJed Brown printf("%s: temperature error %g\n", name, eY.temperature); 34f0b65372SJed Brown PetscFunctionReturn(0); 35f0b65372SJed Brown } 36f0b65372SJed Brown 37f0b65372SJed Brown static PetscErrorCode UnitTests_Newtonian(User user, 38f0b65372SJed Brown NewtonianIdealGasContext gas) { 396dd99bedSLeila Ghaffari 40f0b65372SJed Brown Units units = user->units; 41f0b65372SJed Brown const CeedScalar eps = 1e-6; 42cbe60e31SLeila Ghaffari const CeedScalar kg = units->kilogram, 43cbe60e31SLeila Ghaffari m = units->meter, 44cbe60e31SLeila Ghaffari sec = units->second, 45f0b65372SJed Brown Pascal = units->Pascal; 46f0b65372SJed Brown PetscFunctionBeginUser; 47cbe60e31SLeila Ghaffari const CeedScalar rho = 1.2 * kg / (m*m*m), 48cbe60e31SLeila Ghaffari u = 40 * m/sec; 49f0b65372SJed Brown CeedScalar U[5] = {rho, rho*u, rho *u*1.1, rho *u*1.2, 250e3*Pascal + .5*rho *u*u}; 50f0b65372SJed Brown const CeedScalar x[3] = {.1, .2, .3}; 51f0b65372SJed Brown State s = StateFromU(gas, U, x); 52f0b65372SJed Brown for (int i=0; i<8; i++) { 53f0b65372SJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 54f0b65372SJed Brown if (i < 5) dU[i] = U[i]; 55f0b65372SJed Brown else dx[i-5] = x[i-5]; 56f0b65372SJed Brown State ds = StateFromU_fwd(gas, s, dU, x, dx); 57f0b65372SJed Brown for (int j=0; j<5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 58f0b65372SJed Brown for (int j=0; j<3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 59f0b65372SJed Brown State t = StateFromU(gas, dU, dx); 60f0b65372SJed Brown StatePrimitive dY; 61f0b65372SJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 62f0b65372SJed Brown for (int j=0; j<3; j++) 63f0b65372SJed Brown dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 64f0b65372SJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 65f0b65372SJed Brown char buf[128]; 66f0b65372SJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 67f0b65372SJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 68f0b65372SJed Brown } 69f0b65372SJed Brown PetscFunctionReturn(0); 70f0b65372SJed Brown } 71f0b65372SJed Brown 72b7f03f12SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) { 736dd99bedSLeila Ghaffari 74b7f03f12SJed Brown SetupContext setup_context; 753a8779fbSJames Wright User user = *(User *)ctx; 763a8779fbSJames Wright StabilizationType stab; 773a8779fbSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 783a8779fbSJames Wright PetscBool implicit; 79cbe60e31SLeila Ghaffari PetscBool has_curr_time = PETSC_FALSE, 8034ea8d65SJames Wright use_primitive, unit_tests; 813a8779fbSJames Wright PetscInt ierr; 8215a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 8315a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 843a8779fbSJames Wright 8515a3537eSJed Brown PetscFunctionBeginUser; 86b7f03f12SJed Brown ierr = PetscCalloc1(1, &setup_context); CHKERRQ(ierr); 8715a3537eSJed Brown ierr = PetscCalloc1(1, &newtonian_ig_ctx); CHKERRQ(ierr); 883a8779fbSJames Wright 893a8779fbSJames Wright // ------------------------------------------------------ 903a8779fbSJames Wright // Setup Generic Newtonian IG Problem 913a8779fbSJames Wright // ------------------------------------------------------ 923a8779fbSJames Wright problem->dim = 3; 933a8779fbSJames Wright problem->q_data_size_vol = 10; 94493642f1SJames Wright problem->q_data_size_sur = 10; 95b01ba163SJames Wright problem->jac_data_size_sur = 11; 969785fe93SJed Brown problem->setup_vol.qfunction = Setup; 979785fe93SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 989785fe93SJed Brown problem->setup_sur.qfunction = SetupBoundary; 999785fe93SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 100b7f03f12SJed Brown problem->bc = NULL; 101b7f03f12SJed Brown problem->bc_ctx = setup_context; 1023a8779fbSJames Wright problem->non_zero_time = PETSC_FALSE; 103cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 1043a8779fbSJames Wright 1053a8779fbSJames Wright // ------------------------------------------------------ 1063a8779fbSJames Wright // Create the libCEED context 1073a8779fbSJames Wright // ------------------------------------------------------ 1083a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 1093a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 110bb8a0c61SJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 1113a8779fbSJames Wright CeedScalar lambda = -2./3.; // - 112bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 1133a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 1143a8779fbSJames Wright CeedScalar c_tau = 0.5; // - 115bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 116bb8a0c61SJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 117bb8a0c61SJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 118bb8a0c61SJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 119bb8a0c61SJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 1203a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 1213a8779fbSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 122493642f1SJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 1233a8779fbSJames Wright 1243a8779fbSJames Wright // ------------------------------------------------------ 1253a8779fbSJames Wright // Create the PETSc context 1263a8779fbSJames Wright // ------------------------------------------------------ 127bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 128bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 129bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 1303a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 1313a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 1323a8779fbSJames Wright 1333a8779fbSJames Wright // ------------------------------------------------------ 1343a8779fbSJames Wright // Command line Options 1353a8779fbSJames Wright // ------------------------------------------------------ 1361485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", 1371485969bSJeremy L Thompson NULL); 138cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 139cbe60e31SLeila Ghaffari ierr = PetscOptionsBool("-primitive", "Use primitive variables", 14034ea8d65SJames Wright NULL, use_primitive=PETSC_FALSE, &use_primitive, NULL); CHKERRQ(ierr); 141*d4559bbeSJames Wright // *INDENT-OFF* 14234ea8d65SJames Wright if (use_primitive) { 143cbe60e31SLeila Ghaffari problem->ics.qfunction = ICsNewtonianIG_Prim; 144cbe60e31SLeila Ghaffari problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 145cbe60e31SLeila Ghaffari problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 146cbe60e31SLeila Ghaffari problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 147cbe60e31SLeila Ghaffari problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 148cbe60e31SLeila Ghaffari problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 149*d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 150*d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 151*d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 152*d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 153*d4559bbeSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 154*d4559bbeSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 155*d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 156*d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 157cbe60e31SLeila Ghaffari } else { 158cbe60e31SLeila Ghaffari problem->ics.qfunction = ICsNewtonianIG; 159cbe60e31SLeila Ghaffari problem->ics.qfunction_loc = ICsNewtonianIG_loc; 160cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 161cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 162cbe60e31SLeila Ghaffari problem->apply_vol_ifunction.qfunction = IFunction_Newtonian; 163cbe60e31SLeila Ghaffari problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_loc; 164cbe60e31SLeila Ghaffari problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian; 165cbe60e31SLeila Ghaffari problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_loc; 166*d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 167*d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 168*d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 169*d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 170*d4559bbeSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 171*d4559bbeSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 172*d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 173*d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 174cbe60e31SLeila Ghaffari } 175*d4559bbeSJames Wright // *INDENT-ON* 1761485969bSJeremy L Thompson 1773a8779fbSJames Wright // -- Physics 1783a8779fbSJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 1793a8779fbSJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 1803a8779fbSJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 1813a8779fbSJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 1823a8779fbSJames Wright ierr = PetscOptionsScalar("-lambda", 1833a8779fbSJames Wright "Stokes hypothesis second viscosity coefficient", 1843a8779fbSJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 1853a8779fbSJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 1863a8779fbSJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 1873a8779fbSJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 1883a8779fbSJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 1893a8779fbSJames Wright 190bb8a0c61SJames Wright PetscInt dim = problem->dim; 191bb8a0c61SJames Wright ierr = PetscOptionsRealArray("-g", "Gravitational acceleration", 192bb8a0c61SJames Wright NULL, g, &dim, NULL); CHKERRQ(ierr); 1933a8779fbSJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 1943a8779fbSJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 1953a8779fbSJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 1963a8779fbSJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 1973a8779fbSJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 198bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant", 199bb8a0c61SJames Wright NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr); 200bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", 201bb8a0c61SJames Wright NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr); 202bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", 203bb8a0c61SJames Wright NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr); 204bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", 205bb8a0c61SJames Wright NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr); 206bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", 207bb8a0c61SJames Wright NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr); 2083a8779fbSJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 2093a8779fbSJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 2103a8779fbSJames Wright CHKERRQ(ierr); 211f0b65372SJed Brown ierr = PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", 212f0b65372SJed Brown NULL, unit_tests=PETSC_FALSE, &unit_tests, NULL); 213f0b65372SJed Brown CHKERRQ(ierr); 2143a8779fbSJames Wright 2153a8779fbSJames Wright // -- Units 2163a8779fbSJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 2173a8779fbSJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 2183a8779fbSJames Wright meter = fabs(meter); 2193a8779fbSJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 2203a8779fbSJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 2213a8779fbSJames Wright kilogram = fabs(kilogram); 2223a8779fbSJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 2233a8779fbSJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 2243a8779fbSJames Wright second = fabs(second); 2253a8779fbSJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 2263a8779fbSJames Wright "1 Kelvin in scaled temperature units", 2273a8779fbSJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 2283a8779fbSJames Wright Kelvin = fabs(Kelvin); 2293a8779fbSJames Wright 2303a8779fbSJames Wright // -- Warnings 2313a8779fbSJames Wright if (stab == STAB_SUPG && !implicit) { 2323a8779fbSJames Wright ierr = PetscPrintf(comm, 2333a8779fbSJames Wright "Warning! Use -stab supg only with -implicit\n"); 2343a8779fbSJames Wright CHKERRQ(ierr); 2353a8779fbSJames Wright } 23634ea8d65SJames Wright if (use_primitive && !implicit) { 237cbe60e31SLeila Ghaffari SETERRQ(comm, PETSC_ERR_ARG_NULL, 238cbe60e31SLeila Ghaffari "RHSFunction is not provided for primitive variables (use -primitive only with -implicit)\n"); 239cbe60e31SLeila Ghaffari } 2401485969bSJeremy L Thompson PetscOptionsEnd(); 2413a8779fbSJames Wright 2423a8779fbSJames Wright // ------------------------------------------------------ 2433a8779fbSJames Wright // Set up the PETSc context 2443a8779fbSJames Wright // ------------------------------------------------------ 2453a8779fbSJames Wright // -- Define derived units 2463a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 2473a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 2483a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 2493a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 2503a8779fbSJames Wright 2513a8779fbSJames Wright user->units->meter = meter; 2523a8779fbSJames Wright user->units->kilogram = kilogram; 2533a8779fbSJames Wright user->units->second = second; 2543a8779fbSJames Wright user->units->Kelvin = Kelvin; 2553a8779fbSJames Wright user->units->Pascal = Pascal; 2563a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 2573a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 2583a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 2593a8779fbSJames Wright 2603a8779fbSJames Wright // ------------------------------------------------------ 2613a8779fbSJames Wright // Set up the libCEED context 2623a8779fbSJames Wright // ------------------------------------------------------ 2633a8779fbSJames Wright // -- Scale variables to desired units 2643a8779fbSJames Wright cv *= J_per_kg_K; 2653a8779fbSJames Wright cp *= J_per_kg_K; 2663a8779fbSJames Wright mu *= Pascal * second; 2673a8779fbSJames Wright k *= W_per_m_K; 268493642f1SJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] *= meter; 269493642f1SJames Wright for (PetscInt i=0; i<3; i++) g[i] *= m_per_squared_s; 2703a8779fbSJames Wright problem->dm_scale = meter; 2713a8779fbSJames Wright 2723a8779fbSJames Wright // -- Setup Context 2733a8779fbSJames Wright setup_context->cv = cv; 2743a8779fbSJames Wright setup_context->cp = cp; 2753a8779fbSJames Wright setup_context->lx = domain_size[0]; 2763a8779fbSJames Wright setup_context->ly = domain_size[1]; 2773a8779fbSJames Wright setup_context->lz = domain_size[2]; 2783a8779fbSJames Wright setup_context->time = 0; 279bb8a0c61SJames Wright ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr); 2803a8779fbSJames Wright 2813a8779fbSJames Wright // -- Solver Settings 2823a8779fbSJames Wright user->phys->stab = stab; 2833a8779fbSJames Wright user->phys->implicit = implicit; 28434ea8d65SJames Wright user->phys->use_primitive = use_primitive; 2853a8779fbSJames Wright user->phys->has_curr_time = has_curr_time; 2863a8779fbSJames Wright 2873a8779fbSJames Wright // -- QFunction Context 28815a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 28915a3537eSJed Brown newtonian_ig_ctx->mu = mu; 29015a3537eSJed Brown newtonian_ig_ctx->k = k; 29115a3537eSJed Brown newtonian_ig_ctx->cv = cv; 29215a3537eSJed Brown newtonian_ig_ctx->cp = cp; 29315a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 29415a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 29515a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 29615a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 29715a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 29815a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 29915a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 3008085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 30134ea8d65SJames Wright newtonian_ig_ctx->use_primitive = use_primitive; 30215a3537eSJed Brown ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 3033a8779fbSJames Wright 30415a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 30515a3537eSJed Brown CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, 30615a3537eSJed Brown CEED_USE_POINTER, sizeof(*setup_context), setup_context); 307270bbb13SJeremy L Thompson CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, 308270bbb13SJeremy L Thompson CEED_MEM_HOST, 309270bbb13SJeremy L Thompson FreeContextPetsc); 31015a3537eSJed Brown CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, 31115a3537eSJed Brown "evaluation time", 31215a3537eSJed Brown (char *)&setup_context->time - (char *)setup_context, 1, "Time of evaluation"); 31315a3537eSJed Brown 31415a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 31515a3537eSJed Brown CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, 31615a3537eSJed Brown CEED_USE_POINTER, 31715a3537eSJed Brown sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 31815a3537eSJed Brown CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, 31915a3537eSJed Brown FreeContextPetsc); 32015a3537eSJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", 32115a3537eSJed Brown offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t"); 322f0b65372SJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 323f0b65372SJed Brown offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 324f0b65372SJed Brown "Shift for mass matrix in IJacobian"); 32515a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 32615a3537eSJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 32715a3537eSJed Brown &problem->apply_vol_ifunction.qfunction_context); 328f0b65372SJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 329f0b65372SJed Brown &problem->apply_vol_ijacobian.qfunction_context); 3308085925cSJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 3318085925cSJames Wright &problem->apply_inflow.qfunction_context); 33204b9037bSJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 3330aa41abfSJames Wright &problem->apply_inflow_jacobian.qfunction_context); 3340aa41abfSJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 33504b9037bSJames Wright &problem->apply_outflow.qfunction_context); 33604b9037bSJames Wright CeedQFunctionContextReferenceCopy(newtonian_ig_context, 33704b9037bSJames Wright &problem->apply_outflow_jacobian.qfunction_context); 338f0b65372SJed Brown 339f0b65372SJed Brown if (unit_tests) { 340f0b65372SJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 341f0b65372SJed Brown } 3423a8779fbSJames Wright PetscFunctionReturn(0); 3433a8779fbSJames Wright } 3443a8779fbSJames Wright 345cbe60e31SLeila Ghaffari PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) { 346cbe60e31SLeila Ghaffari 34715a3537eSJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 34815a3537eSJed Brown PetscErrorCode ierr; 34915a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 35015a3537eSJed Brown 3513a8779fbSJames Wright PetscFunctionBeginUser; 35215a3537eSJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 35315a3537eSJed Brown CEED_MEM_HOST, &newtonian_ctx); 35415a3537eSJed Brown ierr = PetscPrintf(comm, 35515a3537eSJed Brown " Problem:\n" 35615a3537eSJed Brown " Problem Name : %s\n" 35715a3537eSJed Brown " Stabilization : %s\n", 35815a3537eSJed Brown app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]); 35915a3537eSJed Brown CHKERRQ(ierr); 36015a3537eSJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 36115a3537eSJed Brown &newtonian_ctx); 3623a8779fbSJames Wright PetscFunctionReturn(0); 3633a8779fbSJames Wright } 364