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 15*e334ad8fSJed Brown // Compute relative error |a - b|/|s| 16*e334ad8fSJed Brown static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, 17*e334ad8fSJed Brown StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 18*e334ad8fSJed Brown PetscReal rtol_velocity, PetscReal rtol_temperature) { 19*e334ad8fSJed Brown PetscFunctionBeginUser; 20*e334ad8fSJed Brown StatePrimitive eY; // relative error 21*e334ad8fSJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 22*e334ad8fSJed Brown PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square( 23*e334ad8fSJed Brown sY.velocity[2])); 24*e334ad8fSJed Brown for (int j=0; j<3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 25*e334ad8fSJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 26*e334ad8fSJed Brown if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, 27*e334ad8fSJed Brown eY.pressure); 28*e334ad8fSJed Brown for (int j=0; j<3; 29*e334ad8fSJed Brown j++) if (fabs(eY.velocity[j]) > rtol_velocity) 30*e334ad8fSJed Brown printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 31*e334ad8fSJed Brown if (fabs(eY.temperature) > rtol_temperature) 32*e334ad8fSJed Brown printf("%s: temperature error %g\n", name, eY.temperature); 33*e334ad8fSJed Brown PetscFunctionReturn(0); 34*e334ad8fSJed Brown } 35*e334ad8fSJed Brown 36*e334ad8fSJed Brown static PetscErrorCode UnitTests_Newtonian(User user, 37*e334ad8fSJed Brown NewtonianIdealGasContext gas) { 38*e334ad8fSJed Brown Units units = user->units; 39*e334ad8fSJed Brown const CeedScalar eps = 1e-6; 40*e334ad8fSJed Brown const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, 41*e334ad8fSJed Brown Pascal = units->Pascal; 42*e334ad8fSJed Brown 43*e334ad8fSJed Brown PetscFunctionBeginUser; 44*e334ad8fSJed Brown const CeedScalar rho = 1.2 * kg / (m*m*m), u = 40 * m/sec; 45*e334ad8fSJed Brown CeedScalar U[5] = {rho, rho*u, rho *u*1.1, rho *u*1.2, 250e3*Pascal + .5*rho *u*u}; 46*e334ad8fSJed Brown const CeedScalar x[3] = {.1, .2, .3}; 47*e334ad8fSJed Brown State s = StateFromU(gas, U, x); 48*e334ad8fSJed Brown for (int i=0; i<8; i++) { 49*e334ad8fSJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 50*e334ad8fSJed Brown if (i < 5) dU[i] = U[i]; 51*e334ad8fSJed Brown else dx[i-5] = x[i-5]; 52*e334ad8fSJed Brown State ds = StateFromU_fwd(gas, s, dU, x, dx); 53*e334ad8fSJed Brown for (int j=0; j<5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 54*e334ad8fSJed Brown for (int j=0; j<3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 55*e334ad8fSJed Brown State t = StateFromU(gas, dU, dx); 56*e334ad8fSJed Brown StatePrimitive dY; 57*e334ad8fSJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 58*e334ad8fSJed Brown for (int j=0; j<3; j++) 59*e334ad8fSJed Brown dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 60*e334ad8fSJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 61*e334ad8fSJed Brown char buf[128]; 62*e334ad8fSJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 63*e334ad8fSJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 64*e334ad8fSJed Brown } 65*e334ad8fSJed Brown PetscFunctionReturn(0); 66*e334ad8fSJed Brown } 67*e334ad8fSJed Brown 68a0add3c9SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) { 69a0add3c9SJed Brown SetupContext setup_context; 7088b783a1SJames Wright User user = *(User *)ctx; 7188b783a1SJames Wright StabilizationType stab; 7288b783a1SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 7388b783a1SJames Wright PetscBool implicit; 74*e334ad8fSJed Brown PetscBool has_curr_time = PETSC_FALSE, unit_tests; 7588b783a1SJames Wright PetscInt ierr; 76841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 77841e4c73SJed Brown CeedQFunctionContext newtonian_ig_context; 7888b783a1SJames Wright 79841e4c73SJed Brown PetscFunctionBeginUser; 80a0add3c9SJed Brown ierr = PetscCalloc1(1, &setup_context); CHKERRQ(ierr); 81841e4c73SJed Brown ierr = PetscCalloc1(1, &newtonian_ig_ctx); CHKERRQ(ierr); 8288b783a1SJames Wright 8388b783a1SJames Wright // ------------------------------------------------------ 8488b783a1SJames Wright // Setup Generic Newtonian IG Problem 8588b783a1SJames Wright // ------------------------------------------------------ 8688b783a1SJames Wright problem->dim = 3; 8788b783a1SJames Wright problem->q_data_size_vol = 10; 8888b783a1SJames Wright problem->q_data_size_sur = 4; 89*e334ad8fSJed Brown problem->jac_data_size_sur = 5; 9091e5af17SJed Brown problem->setup_vol.qfunction = Setup; 9191e5af17SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 9291e5af17SJed Brown problem->ics.qfunction = ICsNewtonianIG; 9391e5af17SJed Brown problem->ics.qfunction_loc = ICsNewtonianIG_loc; 9491e5af17SJed Brown problem->setup_sur.qfunction = SetupBoundary; 9591e5af17SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 965c677226SJed Brown problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 975c677226SJed Brown problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 9891e5af17SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Newtonian; 9991e5af17SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_loc; 100*e334ad8fSJed Brown problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian; 101*e334ad8fSJed Brown problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_loc; 102a0add3c9SJed Brown problem->bc = NULL; 103a0add3c9SJed Brown problem->bc_ctx = setup_context; 10488b783a1SJames Wright problem->non_zero_time = PETSC_FALSE; 10588b783a1SJames Wright problem->print_info = PRINT_DENSITY_CURRENT; 10688b783a1SJames Wright 10788b783a1SJames Wright // ------------------------------------------------------ 10888b783a1SJames Wright // Create the libCEED context 10988b783a1SJames Wright // ------------------------------------------------------ 11088b783a1SJames Wright CeedScalar cv = 717.; // J/(kg K) 11188b783a1SJames Wright CeedScalar cp = 1004.; // J/(kg K) 11288626eedSJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 11388b783a1SJames Wright CeedScalar lambda = -2./3.; // - 11488626eedSJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 11588b783a1SJames Wright CeedScalar k = 0.02638; // W/(m K) 11688b783a1SJames Wright CeedScalar c_tau = 0.5; // - 11788626eedSJames Wright CeedScalar Ctau_t = 1.0; // - 11888626eedSJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 11988626eedSJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 12088626eedSJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 12188626eedSJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 12288b783a1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 12388b783a1SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 12488b783a1SJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 12588b783a1SJames Wright 12688b783a1SJames Wright // ------------------------------------------------------ 12788b783a1SJames Wright // Create the PETSc context 12888b783a1SJames Wright // ------------------------------------------------------ 12988626eedSJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 13088626eedSJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 13188626eedSJames Wright PetscScalar second = 1; // 1 second in scaled time units 13288b783a1SJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 13388b783a1SJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 13488b783a1SJames Wright 13588b783a1SJames Wright // ------------------------------------------------------ 13688b783a1SJames Wright // Command line Options 13788b783a1SJames Wright // ------------------------------------------------------ 13867490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", 13967490bc6SJeremy L Thompson NULL); 14067490bc6SJeremy L Thompson 14188b783a1SJames Wright // -- Physics 14288b783a1SJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 14388b783a1SJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 14488b783a1SJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 14588b783a1SJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 14688b783a1SJames Wright ierr = PetscOptionsScalar("-lambda", 14788b783a1SJames Wright "Stokes hypothesis second viscosity coefficient", 14888b783a1SJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 14988b783a1SJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 15088b783a1SJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 15188b783a1SJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 15288b783a1SJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 15388b783a1SJames Wright 15488626eedSJames Wright PetscInt dim = problem->dim; 15588626eedSJames Wright ierr = PetscOptionsRealArray("-g", "Gravitational acceleration", 15688626eedSJames Wright NULL, g, &dim, NULL); CHKERRQ(ierr); 15788b783a1SJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 15888b783a1SJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 15988b783a1SJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 16088b783a1SJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 16188b783a1SJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 16288626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant", 16388626eedSJames Wright NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr); 16488626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", 16588626eedSJames Wright NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr); 16688626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", 16788626eedSJames Wright NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr); 16888626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", 16988626eedSJames Wright NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr); 17088626eedSJames Wright ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", 17188626eedSJames Wright NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr); 17288b783a1SJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 17388b783a1SJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 17488b783a1SJames Wright CHKERRQ(ierr); 175*e334ad8fSJed Brown ierr = PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", 176*e334ad8fSJed Brown NULL, unit_tests=PETSC_FALSE, &unit_tests, NULL); 177*e334ad8fSJed Brown CHKERRQ(ierr); 17888b783a1SJames Wright 17988b783a1SJames Wright // -- Units 18088b783a1SJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 18188b783a1SJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 18288b783a1SJames Wright meter = fabs(meter); 18388b783a1SJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 18488b783a1SJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 18588b783a1SJames Wright kilogram = fabs(kilogram); 18688b783a1SJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 18788b783a1SJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 18888b783a1SJames Wright second = fabs(second); 18988b783a1SJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 19088b783a1SJames Wright "1 Kelvin in scaled temperature units", 19188b783a1SJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 19288b783a1SJames Wright Kelvin = fabs(Kelvin); 19388b783a1SJames Wright 19488b783a1SJames Wright // -- Warnings 19588b783a1SJames Wright if (stab == STAB_SUPG && !implicit) { 19688b783a1SJames Wright ierr = PetscPrintf(comm, 19788b783a1SJames Wright "Warning! Use -stab supg only with -implicit\n"); 19888b783a1SJames Wright CHKERRQ(ierr); 19988b783a1SJames Wright } 20067490bc6SJeremy L Thompson PetscOptionsEnd(); 20188b783a1SJames Wright 20288b783a1SJames Wright // ------------------------------------------------------ 20388b783a1SJames Wright // Set up the PETSc context 20488b783a1SJames Wright // ------------------------------------------------------ 20588b783a1SJames Wright // -- Define derived units 20688b783a1SJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 20788b783a1SJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 20888b783a1SJames Wright m_per_squared_s = meter / PetscSqr(second); 20988b783a1SJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 21088b783a1SJames Wright 21188b783a1SJames Wright user->units->meter = meter; 21288b783a1SJames Wright user->units->kilogram = kilogram; 21388b783a1SJames Wright user->units->second = second; 21488b783a1SJames Wright user->units->Kelvin = Kelvin; 21588b783a1SJames Wright user->units->Pascal = Pascal; 21688b783a1SJames Wright user->units->J_per_kg_K = J_per_kg_K; 21788b783a1SJames Wright user->units->m_per_squared_s = m_per_squared_s; 21888b783a1SJames Wright user->units->W_per_m_K = W_per_m_K; 21988b783a1SJames Wright 22088b783a1SJames Wright // ------------------------------------------------------ 22188b783a1SJames Wright // Set up the libCEED context 22288b783a1SJames Wright // ------------------------------------------------------ 22388b783a1SJames Wright // -- Scale variables to desired units 22488b783a1SJames Wright cv *= J_per_kg_K; 22588b783a1SJames Wright cp *= J_per_kg_K; 22688b783a1SJames Wright mu *= Pascal * second; 22788b783a1SJames Wright k *= W_per_m_K; 22888b783a1SJames Wright for (int i=0; i<3; i++) domain_size[i] *= meter; 22988626eedSJames Wright for (int i=0; i<3; i++) g[i] *= m_per_squared_s; 23088b783a1SJames Wright problem->dm_scale = meter; 23188b783a1SJames Wright 23288b783a1SJames Wright // -- Setup Context 23388b783a1SJames Wright setup_context->cv = cv; 23488b783a1SJames Wright setup_context->cp = cp; 23588b783a1SJames Wright setup_context->lx = domain_size[0]; 23688b783a1SJames Wright setup_context->ly = domain_size[1]; 23788b783a1SJames Wright setup_context->lz = domain_size[2]; 23888b783a1SJames Wright setup_context->time = 0; 23988626eedSJames Wright ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr); 24088b783a1SJames Wright 24188b783a1SJames Wright // -- Solver Settings 24288b783a1SJames Wright user->phys->stab = stab; 24388b783a1SJames Wright user->phys->implicit = implicit; 24488b783a1SJames Wright user->phys->has_curr_time = has_curr_time; 24588b783a1SJames Wright 24688b783a1SJames Wright // -- QFunction Context 247841e4c73SJed Brown newtonian_ig_ctx->lambda = lambda; 248841e4c73SJed Brown newtonian_ig_ctx->mu = mu; 249841e4c73SJed Brown newtonian_ig_ctx->k = k; 250841e4c73SJed Brown newtonian_ig_ctx->cv = cv; 251841e4c73SJed Brown newtonian_ig_ctx->cp = cp; 252841e4c73SJed Brown newtonian_ig_ctx->c_tau = c_tau; 253841e4c73SJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 254841e4c73SJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 255841e4c73SJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 256841e4c73SJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 257841e4c73SJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 258841e4c73SJed Brown newtonian_ig_ctx->stabilization = stab; 259841e4c73SJed Brown ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 26088b783a1SJames Wright 261841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 262841e4c73SJed Brown CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, 263841e4c73SJed Brown CEED_USE_POINTER, sizeof(*setup_context), setup_context); 264841e4c73SJed Brown CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, 265841e4c73SJed Brown "evaluation time", 266841e4c73SJed Brown (char *)&setup_context->time - (char *)setup_context, 1, "Time of evaluation"); 267841e4c73SJed Brown 268841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 269841e4c73SJed Brown CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, 270841e4c73SJed Brown CEED_USE_POINTER, 271841e4c73SJed Brown sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 272841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, 273841e4c73SJed Brown FreeContextPetsc); 274841e4c73SJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", 275841e4c73SJed Brown offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t"); 276*e334ad8fSJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 277*e334ad8fSJed Brown offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 278*e334ad8fSJed Brown "Shift for mass matrix in IJacobian"); 279841e4c73SJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 280841e4c73SJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 281841e4c73SJed Brown &problem->apply_vol_ifunction.qfunction_context); 282*e334ad8fSJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 283*e334ad8fSJed Brown &problem->apply_vol_ijacobian.qfunction_context); 284*e334ad8fSJed Brown 285*e334ad8fSJed Brown if (unit_tests) { 286*e334ad8fSJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 287*e334ad8fSJed Brown } 28888b783a1SJames Wright PetscFunctionReturn(0); 28988b783a1SJames Wright } 29088b783a1SJames Wright 291841e4c73SJed Brown PetscErrorCode PRINT_DENSITY_CURRENT(ProblemData *problem, 292841e4c73SJed Brown AppCtx app_ctx) { 293841e4c73SJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 294841e4c73SJed Brown PetscErrorCode ierr; 295841e4c73SJed Brown NewtonianIdealGasContext newtonian_ctx; 296841e4c73SJed Brown 29788b783a1SJames Wright PetscFunctionBeginUser; 298841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 299841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ctx); 300841e4c73SJed Brown ierr = PetscPrintf(comm, 301841e4c73SJed Brown " Problem:\n" 302841e4c73SJed Brown " Problem Name : %s\n" 303841e4c73SJed Brown " Stabilization : %s\n", 304841e4c73SJed Brown app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]); 305841e4c73SJed Brown CHKERRQ(ierr); 306841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 307841e4c73SJed Brown &newtonian_ctx); 30888b783a1SJames Wright PetscFunctionReturn(0); 30988b783a1SJames Wright } 310