Lines Matching +full:- +full:o

1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
20 MPI_Comm comm = honee->comm; in PRINT_NEWTONIAN()
21 Ceed ceed = honee->ceed; in PRINT_NEWTONIAN()
25 …PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newt… in PRINT_NEWTONIAN()
30 app_ctx->problem_name, StabilizationTypes[newt_ctx->stabilization])); in PRINT_NEWTONIAN()
31 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newt_ctx)); in PRINT_NEWTONIAN()
39 Ceed ceed = honee->ceed; in CreateKSPMassOperator_NewtonianStabilized()
54 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); in CreateKSPMassOperator_NewtonianStabilized()
79 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); in CreateKSPMassOperator_NewtonianStabilized()
101 Ceed ceed = honee->ceed; in DivDiffFluxProjectionCreateRHS_Direct_NS()
102 NodalProjectionData projection = diff_flux_proj->projection; in DivDiffFluxProjectionCreateRHS_Direct_NS()
109 // -- Get Pre-requisite things in DivDiffFluxProjectionCreateRHS_Direct_NS()
110 PetscCall(DMGetDimension(projection->dm, &dim)); in DivDiffFluxProjectionCreateRHS_Direct_NS()
111 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q)); in DivDiffFluxProjectionCreateRHS_Direct_NS()
114 …CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op; in DivDiffFluxProjectionCreateRHS_Direct_NS()
130 …PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE,… in DivDiffFluxProjectionCreateRHS_Direct_NS()
131 …PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &ba… in DivDiffFluxProjectionCreateRHS_Direct_NS()
132 …PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &… in DivDiffFluxProjectionCreateRHS_Direct_NS()
133 switch (honee->phys->state_var) { in DivDiffFluxProjectionCreateRHS_Direct_NS()
152 …PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_co… in DivDiffFluxProjectionCreateRHS_Direct_NS()
178 // -- Build RHS operator in DivDiffFluxProjectionCreateRHS_Direct_NS()
179 switch (honee->phys->state_var) { in DivDiffFluxProjectionCreateRHS_Direct_NS()
194 PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); in DivDiffFluxProjectionCreateRHS_Direct_NS()
199 …PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_… in DivDiffFluxProjectionCreateRHS_Direct_NS()
201 PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); in DivDiffFluxProjectionCreateRHS_Direct_NS()
202 …PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &fa… in DivDiffFluxProjectionCreateRHS_Direct_NS()
210 …PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); in DivDiffFluxProjectionCreateRHS_Direct_NS()
211 PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); in DivDiffFluxProjectionCreateRHS_Direct_NS()
214 …PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_… in DivDiffFluxProjectionCreateRHS_Direct_NS()
215 for (PetscInt o = 0; o < num_orientations_values; o++) { in DivDiffFluxProjectionCreateRHS_Direct_NS() local
221 …PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_f… in DivDiffFluxProjectionCreateRHS_Direct_NS()
223 …PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, he… in DivDiffFluxProjectionCreateRHS_Direct_NS()
224 …PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, or… in DivDiffFluxProjectionCreateRHS_Direct_NS()
225 …PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientatio… in DivDiffFluxProjectionCreateRHS_Direct_NS()
227 …PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_… in DivDiffFluxProjectionCreateRHS_Direct_NS()
228 …PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, &elem_res… in DivDiffFluxProjectionCreateRHS_Direct_NS()
265 Ceed ceed = honee->ceed; in DivDiffFluxProjectionCreateRHS_Indirect_NS()
266 NodalProjectionData projection = diff_flux_proj->projection; in DivDiffFluxProjectionCreateRHS_Indirect_NS()
277 PetscCall(DMGetDimension(projection->dm, &dim)); in DivDiffFluxProjectionCreateRHS_Indirect_NS()
278 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q)); in DivDiffFluxProjectionCreateRHS_Indirect_NS()
281 …CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op; in DivDiffFluxProjectionCreateRHS_Indirect_NS()
287 …PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE,… in DivDiffFluxProjectionCreateRHS_Indirect_NS()
288 …PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &ba… in DivDiffFluxProjectionCreateRHS_Indirect_NS()
289 …PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_V… in DivDiffFluxProjectionCreateRHS_Indirect_NS()
290 …PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, heig… in DivDiffFluxProjectionCreateRHS_Indirect_NS()
291 …PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &… in DivDiffFluxProjectionCreateRHS_Indirect_NS()
293 switch (honee->phys->state_var) { in DivDiffFluxProjectionCreateRHS_Indirect_NS()
309 …PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_I… in DivDiffFluxProjectionCreateRHS_Indirect_NS()
333 Honee honee = honee_bc->honee; in BoundaryIntegralBCSetup_CreateIFunctionQF()
335 switch (honee->phys->state_var) { in BoundaryIntegralBCSetup_CreateIFunctionQF()
337 …eIFunctionQF(bc_def, BoundaryIntegral_Conserv, BoundaryIntegral_Conserv_loc, honee_bc->qfctx, qf)); in BoundaryIntegralBCSetup_CreateIFunctionQF()
340 …CCreateIFunctionQF(bc_def, BoundaryIntegral_Prim, BoundaryIntegral_Prim_loc, honee_bc->qfctx, qf)); in BoundaryIntegralBCSetup_CreateIFunctionQF()
343 …eIFunctionQF(bc_def, BoundaryIntegral_Entropy, BoundaryIntegral_Entropy_loc, honee_bc->qfctx, qf)); in BoundaryIntegralBCSetup_CreateIFunctionQF()
354 Honee honee = honee_bc->honee; in BoundaryIntegralBCSetup_CreateIJacobianQF()
356 switch (honee->phys->state_var) { in BoundaryIntegralBCSetup_CreateIJacobianQF()
358 …f, BoundaryIntegral_Jacobian_Conserv, BoundaryIntegral_Jacobian_Conserv_loc, honee_bc->qfctx, qf)); in BoundaryIntegralBCSetup_CreateIJacobianQF()
361 …(bc_def, BoundaryIntegral_Jacobian_Prim, BoundaryIntegral_Jacobian_Prim_loc, honee_bc->qfctx, qf)); in BoundaryIntegralBCSetup_CreateIJacobianQF()
364 …f, BoundaryIntegral_Jacobian_Entropy, BoundaryIntegral_Jacobian_Entropy_loc, honee_bc->qfctx, qf)); in BoundaryIntegralBCSetup_CreateIJacobianQF()
373 CeedInt degree = honee->app_ctx->degree; in NS_NEWTONIAN_IG()
376 MPI_Comm comm = honee->comm; in NS_NEWTONIAN_IG()
377 Ceed ceed = honee->ceed; in NS_NEWTONIAN_IG()
386 CeedScalar idl_decay_time = -1; in NS_NEWTONIAN_IG()
393 .lambda = -2. / 3., in NS_NEWTONIAN_IG()
394 .mu = 1.8e-5, in NS_NEWTONIAN_IG()
400 .Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1], in NS_NEWTONIAN_IG()
413 …PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)… in NS_NEWTONIAN_IG()
417 …PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, newtonian_ig_ctx->ga… in NS_NEWTONIAN_IG()
418 …etscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, newtonian_ig_ctx->g… in NS_NEWTONIAN_IG()
419 …PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, ne… in NS_NEWTONIAN_IG()
420 &newtonian_ig_ctx->gas.lambda, NULL)); in NS_NEWTONIAN_IG()
421 …etscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, newtonian_ig_ctx->… in NS_NEWTONIAN_IG()
422 …PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, newtonian_ig_ctx->gas.k, &newtoni… in NS_NEWTONIAN_IG()
426 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); in NS_NEWTONIAN_IG()
427 …PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, newtonian_i… in NS_NEWTONIAN_IG()
431 …PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(… in NS_NEWTONIAN_IG()
432 …PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, newtonian_ig_ctx->tau… in NS_NEWTONIAN_IG()
433 &newtonian_ig_ctx->tau_coeffs.Ctau_t, NULL)); in NS_NEWTONIAN_IG()
434 …PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, newtonian_ig_ctx->… in NS_NEWTONIAN_IG()
435 &newtonian_ig_ctx->tau_coeffs.Ctau_v, NULL)); in NS_NEWTONIAN_IG()
436 …PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, newtonian_ig_ct… in NS_NEWTONIAN_IG()
437 &newtonian_ig_ctx->tau_coeffs.Ctau_C, NULL)); in NS_NEWTONIAN_IG()
438 …PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, newtonian_ig_ctx- in NS_NEWTONIAN_IG()
439 &newtonian_ig_ctx->tau_coeffs.Ctau_M, NULL)); in NS_NEWTONIAN_IG()
440 …PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, newtonian_ig_ctx->t… in NS_NEWTONIAN_IG()
441 &newtonian_ig_ctx->tau_coeffs.Ctau_E, NULL)); in NS_NEWTONIAN_IG()
444 …PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.… in NS_NEWTONIAN_IG()
445 …PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, refer… in NS_NEWTONIAN_IG()
446 …PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, refe… in NS_NEWTONIAN_IG()
448 …PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = P… in NS_NEWTONIAN_IG()
449 …PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests =… in NS_NEWTONIAN_IG()
451 …"RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit… in NS_NEWTONIAN_IG()
455 …PetscBool idl_enable = (PetscBool)newtonian_ig_ctx->idl_enable; // Need PetscBool variable to rea… in NS_NEWTONIAN_IG()
456 …PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance… in NS_NEWTONIAN_IG()
460 newtonian_ig_ctx->idl_enable = idl_enable; in NS_NEWTONIAN_IG()
462 …cCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, newtonian_ig_ctx->… in NS_NEWTONIAN_IG()
464 …PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, new… in NS_NEWTONIAN_IG()
465 &newtonian_ig_ctx->idl_length, NULL)); in NS_NEWTONIAN_IG()
466 newtonian_ig_ctx->idl_pressure = reference.pressure; in NS_NEWTONIAN_IG()
467 …PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-refere… in NS_NEWTONIAN_IG()
468 … newtonian_ig_ctx->idl_pressure, &newtonian_ig_ctx->idl_pressure, NULL)); in NS_NEWTONIAN_IG()
472 // ------------------------------------------------------ in NS_NEWTONIAN_IG()
474 // ------------------------------------------------------ in NS_NEWTONIAN_IG()
475 // -- Scale variables to desired units in NS_NEWTONIAN_IG()
476 Units units = honee->units; in NS_NEWTONIAN_IG()
477 newtonian_ig_ctx->gas.cv *= units->J_per_kg_K; in NS_NEWTONIAN_IG()
478 newtonian_ig_ctx->gas.cp *= units->J_per_kg_K; in NS_NEWTONIAN_IG()
479 newtonian_ig_ctx->gas.mu *= units->Pascal * units->second; in NS_NEWTONIAN_IG()
480 newtonian_ig_ctx->gas.k *= units->W_per_m_K; in NS_NEWTONIAN_IG()
481 for (PetscInt i = 0; i < 3; i++) newtonian_ig_ctx->g[i] *= units->m_per_squared_s; in NS_NEWTONIAN_IG()
482 reference.pressure *= units->Pascal; in NS_NEWTONIAN_IG()
483 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= units->meter / units->second; in NS_NEWTONIAN_IG()
484 reference.temperature *= units->Kelvin; in NS_NEWTONIAN_IG()
488 for (PetscInt i = 0; i < 3; i++) domain_size[i] = (domain_max[i] - domain_min[i]); in NS_NEWTONIAN_IG()
490 // -- Solver Settings in NS_NEWTONIAN_IG()
491 honee->phys->implicit = implicit; in NS_NEWTONIAN_IG()
492 honee->phys->state_var = state_var; in NS_NEWTONIAN_IG()
494 // -- QFunction Context in NS_NEWTONIAN_IG()
495 newtonian_ig_ctx->stabilization = stab; in NS_NEWTONIAN_IG()
496 newtonian_ig_ctx->is_implicit = implicit; in NS_NEWTONIAN_IG()
497 newtonian_ig_ctx->state_var = state_var; in NS_NEWTONIAN_IG()
498 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * units->second); in NS_NEWTONIAN_IG()
499 newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; in NS_NEWTONIAN_IG()
501 // -- Setup Context in NS_NEWTONIAN_IG()
514 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); in NS_NEWTONIAN_IG()
520 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx)); in NS_NEWTONIAN_IG()
532 problem->num_comps_jac_data = 14; in NS_NEWTONIAN_IG()
533 if (newtonian_ig_ctx->idl_enable) problem->num_comps_jac_data += 1; in NS_NEWTONIAN_IG()
534 problem->compute_exact_solution_error = PETSC_FALSE; in NS_NEWTONIAN_IG()
535 problem->print_info = PRINT_NEWTONIAN; in NS_NEWTONIAN_IG()
536 problem->num_components = 5; in NS_NEWTONIAN_IG()
537 PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); in NS_NEWTONIAN_IG()
545 …problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Conserv, .qf_loc = ICsN… in NS_NEWTONIAN_IG()
546 …problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHSFunction_Newtonian, .qf_loc = RHSFu… in NS_NEWTONIAN_IG()
547 …problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Conserv, .qf_loc =… in NS_NEWTONIAN_IG()
548 …problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Conserv, .qf_loc =… in NS_NEWTONIAN_IG()
549 …; i < 5; i++) PetscCall(PetscStrallocpy(conserv_component_names[i], &problem->component_names[i])); in NS_NEWTONIAN_IG()
552 …problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Prim, .qf_loc = ICsNewt… in NS_NEWTONIAN_IG()
553 …problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Prim, .qf_loc = IF… in NS_NEWTONIAN_IG()
554 …problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Prim, .qf_loc = IJ… in NS_NEWTONIAN_IG()
555 …for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(prim_component_names[i], &problem->comp… in NS_NEWTONIAN_IG()
558 …problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsNewtonianIG_Entropy, .qf_loc = ICsN… in NS_NEWTONIAN_IG()
559 …problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Newtonian_Entropy, .qf_loc =… in NS_NEWTONIAN_IG()
560 …problem->apply_vol_ijacobian = (HoneeQFSpec){.qf_func_ptr = IJacobian_Newtonian_Entropy, .qf_loc =… in NS_NEWTONIAN_IG()
561 …; i < 5; i++) PetscCall(PetscStrallocpy(entropy_component_names[i], &problem->component_names[i])); in NS_NEWTONIAN_IG()
565 problem->ics.qfctx = ics_qfctx; in NS_NEWTONIAN_IG()
566 problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; in NS_NEWTONIAN_IG()
567 …PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifun… in NS_NEWTONIAN_IG()
568 …PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijac… in NS_NEWTONIAN_IG()
570 …if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_Newtonia… in NS_NEWTONIAN_IG()
572 …PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 4, &honee->diff_… in NS_NEWTONIAN_IG()
573 if (honee->diff_flux_proj) { in NS_NEWTONIAN_IG()
574 DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; in NS_NEWTONIAN_IG()
575 NodalProjectionData projection = diff_flux_proj->projection; in NS_NEWTONIAN_IG()
578 diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_NS; in NS_NEWTONIAN_IG()
579 diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS; in NS_NEWTONIAN_IG()
580 PetscCall(DMGetLocalSection(projection->dm, &section)); in NS_NEWTONIAN_IG()
581 switch (honee->diff_flux_proj->method) { in NS_NEWTONIAN_IG()
605 …SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with … in NS_NEWTONIAN_IG()
606 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); in NS_NEWTONIAN_IG()
611 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { in NS_NEWTONIAN_IG()
612 BCDefinition bc_def = problem->bc_defs[b]; in NS_NEWTONIAN_IG()
626 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &honee_bc->qfctx)); in NS_NEWTONIAN_IG()
627 honee_bc->honee = honee; in NS_NEWTONIAN_IG()
628 honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0; in NS_NEWTONIAN_IG()
636 if (unit_tests) PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx->gas)); in NS_NEWTONIAN_IG()
640 //------------------------------------
642 //------------------------------------
650 relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1); in CheckQWithTolerance()
651 relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1); in CheckQWithTolerance()
656 relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor; in CheckQWithTolerance()
673 // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
689 snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial); in TestState()
697 CeedScalar eps = 4e-7; // Finite difference step in TestState_fwd()
723 for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps; in TestState_fwd()
726 …snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, … in TestState_fwd()
734 Units units = honee->units; in UnitTests_Newtonian()
735 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin; in UnitTests_Newtonian()
741 const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas.cv * T; in UnitTests_Newtonian()
748 const CeedScalar entropy = log(P) - gamma * log(rho); in UnitTests_Newtonian()
752 …const CeedScalar V0[5] = {(gamma - entropy) / (gamma - 1) - rho_div_p * (e_kinetic), rho_div_… in UnitTests_Newtonian()
753 -rho_div_p}; in UnitTests_Newtonian()
763 rtol = 5e-6; in UnitTests_Newtonian()