PWR-HYBRID-3/code/scripts/sim/validate_pj.jl
Dane Sabo c5133401e0 Session work scratch: scram X_exit refactor, hot-standby SOS, fat scram tubes, model cheatsheet, journal entry
Multi-session work bundle on a draft branch.  Splits into a clean
sequence of commits later; pushed here so it isn't lost on a reboot.

Reach work
- code/scripts/reach/reach_scram_pj.jl: shutdown_margin halfspace
  X_exit (replaces "n <= 1e-4 AND T_f bound" framing); per-step
  envelope extraction added.
- code/scripts/reach/reach_scram_pj_fat.jl: per-step envelope
  extraction added; shutdown_margin discharge logic mirrored from the
  tight scram script.  3 probes (10/30/60s) all discharge from the
  fat union polytope.
- code/scripts/reach/reach_scram_full_fat.jl (NEW): full nonlinear
  PKE scram reach with fat entry.  Hits the stiffness wall at
  ~1.5 s plant time as expected; saves NaN-tolerant per-step
  envelopes.  Demonstrates concretely why PJ is the right tool for
  the longer-horizon proof.
- code/scripts/reach/reach_heatup_pj.jl: T_REF_START_C constant
  (entry-conditioned ramp) replaces T_STANDBY-init that was making
  the FL controller command cooling at t=0.  Per-step extraction
  already in place.
- code/configs/heatup/tight.toml: bumped maxsteps; probe horizon
  parameterized.

Hot-standby SOS barrier
- code/scripts/barrier/barrier_sos_2d_shutdown.jl (NEW): mirrors the
  operation SOS machinery on the hot-standby thermal projection.
  Includes the eps-slack pattern (so feasibility doesn't silently
  collapse to B == 0).
- code/scripts/barrier/barrier_sos_2d.jl: refactored to use the same
  helper.
- code/src/sos_barrier.jl (NEW): solve_sos_barrier_2d helper module
  factoring out the SOS construction; eps-slack with eps_cap=1.0 to
  avoid unbounded primal.

Library
- code/src/pke_states.jl (NEW): single source of truth for canonical
  initial-condition vectors per DRC mode (op, shutdown, heatup) keyed
  off plant + predicates.
- code/scripts/sim/{main_mode_sweep,validate_pj}.jl, code/CLAUDE.md:
  migrated to pke_states.

Predicates + invariants
- reachability/predicates.json: new shutdown_margin predicate (1%
  dk/k tech-spec floor, expressed as alpha_f*T_f + alpha_c*T_c
  halfspace).  Used as scram X_exit.

Plot script
- code/scripts/plot/plot_reach_tubes.jl: plot_tubes_scram_pj() with
  variant=:fat|:tight knob; plot_tubes_scram_full() for full-PKE
  3-panel (T_c, T_f, rho); plot_tubes_heatup_pj() reads results/
  not reachability/.

Journal + memory
- journal/entries/2026-04-27-shutdown-sos-and-scram-X_exit.tex (NEW):
  long-form entry on the SOS hot-standby barrier and the scram X_exit
  refactor.
- journal/journal.tex: input chain updated.
- claude_memory/ — three new session notes:
  * 2026-04-27-scram-X_exit-shutdown-margin.md
  * 2026-04-28-DICE-2026-conference-intel.md (people, sessions,
    strategic notes for the May 12 talk)
  * 2026-04-28-path1-sos-pj-sketch.md (sketch of nonlinear-SOS via
    polynomial multiply-through; saved for an overnight session)

Docs
- docs/model_cheatsheet.md (NEW): one-page reference of state vector,
  dynamics, constants, modes, predicates, sanity numbers — the talk
  prep cheatsheet Dane asked for.
- docs/figures/reach_*_tubes.png: regenerated with the new mat data.
- presentations/prelim-presentation/outline.md: revised arc per the
  April-28 review pass (cuts: Lyapunov-fails standalone slide,
  operation-tube standalone slide, SOS standalone; adds: scopes-of-
  control framing, scram on the headline result slide).
- app/predicate_explorer.jl: minor.

Hacker-Split: end-of-session scratch bundle
2026-05-02 23:02:50 -04:00

114 lines
4.5 KiB
Julia

#!/usr/bin/env julia
#
# validate_pj.jl — quantify the prompt-jump approximation error.
#
# Two parallel sims of the same heatup scenario:
# (i) full 10-state PKE+T/H (pke_th_rhs!)
# (ii) reduced 9-state prompt-jump model (pke_th_rhs_pj!)
#
# After the prompt transient (~few hundred microseconds) the two
# trajectories should agree closely on T_c, T_f, T_cold, and on the
# reconstructed n vs the dynamic n. Quantify the error explicitly.
using Pkg
Pkg.activate(joinpath(@__DIR__, "..", ".."))
using Printf
using LinearAlgebra
using OrdinaryDiffEq
using Plots
include(joinpath(@__DIR__, "..", "..", "src", "pke_params.jl"))
include(joinpath(@__DIR__, "..", "..", "src", "pke_states.jl"))
include(joinpath(@__DIR__, "..", "..", "src", "pke_th_rhs.jl"))
include(joinpath(@__DIR__, "..", "..", "src", "pke_th_rhs_pj.jl"))
include(joinpath(@__DIR__, "..", "..", "controllers", "controllers.jl"))
plant = pke_params()
states = pke_states(plant) # falls back to default offset; no predicates needed for validation
T_standby = states.T_standby
# Heatup scenario — same as sim_heatup.jl.
ref_heat = (; T_start=T_standby, T_target=plant.T_c0, ramp_rate=28/3600)
Q_sg = t -> 0.0
# Initial conditions — both at n=1e-3 (the heatup canonical IC), same temperatures.
x0_full = states.heatup # 10-state
x0_pj = states.heatup[2:end] # drop n; PJ reconstructs it
tspan = (0.0, 3000.0)
# --- Full-state sim ---
function rhs_full!(dx, x, p, t)
u = ctrl_heatup(t, x, plant, ref_heat)
pke_th_rhs!(dx, x, t, plant, Q_sg, u)
end
prob_full = ODEProblem(rhs_full!, x0_full, tspan)
sol_full = solve(prob_full, Rodas5(); reltol=1e-8, abstol=1e-10)
# --- Prompt-jump sim ---
# ctrl_heatup expects 10-state x with T_f at index 8, T_c at 9.
# We pass an adapter that maps 9-state to the controller's expected layout.
function ctrl_heatup_pj(t, x_pj, plant_arg, ref_arg)
# Construct a fake 10-state with n placeholder; controller doesn't read n.
x10 = [0.0; x_pj[1:6]; x_pj[7]; x_pj[8]; x_pj[9]]
return ctrl_heatup(t, x10, plant_arg, ref_arg)
end
function rhs_pj!(dx, x, p, t)
u = ctrl_heatup_pj(t, x, plant, ref_heat)
pke_th_rhs_pj!(dx, x, t, plant, Q_sg, u)
end
prob_pj = ODEProblem(rhs_pj!, x0_pj, tspan)
sol_pj = solve(prob_pj, Rodas5(); reltol=1e-8, abstol=1e-10)
# --- Compare at sampled times ---
function get_n_full(sol, t) sol(t)[1] end
function get_T_c_full(sol, t) sol(t)[9] end
function get_T_f_full(sol, t) sol(t)[8] end
function get_T_cold_full(sol, t) sol(t)[10] end
function get_T_c_pj(sol, t) sol(t)[8] end
function get_T_f_pj(sol, t) sol(t)[7] end
function get_T_cold_pj(sol, t) sol(t)[9] end
function get_n_pj(sol, t)
x = sol(t)
u = ctrl_heatup_pj(t, x, plant, ref_heat)
pj_reconstruct_n(x, plant, u)
end
println("\n=== PJ vs full-state, heatup scenario ===")
println(" t [s] n_full n_pj |Δn|/n_full T_c err T_f err T_cold err")
for t_check in (1.0, 5.0, 10.0, 60.0, 300.0, 1200.0, 3000.0)
n_f = get_n_full(sol_full, t_check)
n_p = get_n_pj(sol_pj, t_check)
Tc_err = abs(get_T_c_full(sol_full, t_check) - get_T_c_pj(sol_pj, t_check))
Tf_err = abs(get_T_f_full(sol_full, t_check) - get_T_f_pj(sol_pj, t_check))
Tcold_err = abs(get_T_cold_full(sol_full, t_check) - get_T_cold_pj(sol_pj, t_check))
@printf " %6.1f %.3e %.3e %8.2e %.3e %.3e %.3e\n" t_check n_f n_p abs(n_f-n_p)/abs(n_f) Tc_err Tf_err Tcold_err
end
# --- Plot trajectory overlay ---
figdir = joinpath(@__DIR__, "..", "..", "..", "docs", "figures")
isdir(figdir) || mkpath(figdir)
CtoF(T) = T*9/5 + 32
t_grid = range(0, 3000, length=600)
n_full_arr = [get_n_full(sol_full, t) for t in t_grid]
n_pj_arr = [get_n_pj(sol_pj, t) for t in t_grid]
Tc_full_arr = [get_T_c_full(sol_full, t) for t in t_grid]
Tc_pj_arr = [get_T_c_pj(sol_pj, t) for t in t_grid]
p_n = plot(t_grid ./ 60, n_full_arr, lw=2, color=:blue, label="full-state",
xlabel="t [min]", ylabel="n", title="Power: full vs PJ")
plot!(p_n, t_grid ./ 60, n_pj_arr, lw=1.5, ls=:dash, color=:red, label="prompt-jump")
p_Tc = plot(t_grid ./ 60, CtoF.(Tc_full_arr), lw=2, color=:blue, label="full-state",
xlabel="t [min]", ylabel="T_avg [F]", title="T_avg: full vs PJ")
plot!(p_Tc, t_grid ./ 60, CtoF.(Tc_pj_arr), lw=1.5, ls=:dash, color=:red, label="prompt-jump")
fig = plot(p_n, p_Tc, layout=(1,2), size=(1200, 450),
plot_title="PJ vs full-state, heatup scenario (3000 s)")
savefig(fig, joinpath(figdir, "validate_pj_heatup.png"))
println("\nFigure: $figdir/validate_pj_heatup.png")