1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
use std::{
collections::{BTreeSet, HashMap},
time::Instant,
};
use super::diagnostics::{solve_timing_enabled, SparsePhaseTimings};
use super::sparse_order::{build_csc, genmmd_order, greedy_mdo, symbolic_fill};
/// Sparse symmetric positive-definite solver (§3.6).
///
/// Implements the three-phase algorithm described in §3.6:
/// 1. MMD reordering (done once in `new`)
/// 2. Symbolic factorisation (done once in `new`)
/// 3. Numerical factorisation + solve (`factorize_solve`, called every Newton step)
///
/// Matrix values are assembled externally: call `clear` to zero all entries,
/// set `aii[row[ji]]` and `aij[link_aij_pos[k]]` for each junction and link,
/// then call `factorize_solve` which overwrites `aii`/`aij` with L and stores
/// the solution in `f`.
pub struct SparseSolver {
/// Number of junctions (system size).
pub(crate) n: usize,
/// Total below-diagonal non-zeros in L.
#[allow(dead_code)]
pub(crate) n_coeff: usize,
/// Inverse permutation: `row[i]` = elimination step for original junction i.
pub(crate) row: Vec<usize>,
/// Column pointer: `xlnz[k]` = start of column k in `nzsub` (0-based).
pub(crate) xlnz: Vec<usize>,
/// Row indices of non-zeros (permuted ordering, sorted per column).
pub(crate) nzsub: Vec<usize>,
/// Off-diagonal values, size `n_coeff`. Overwritten in-place during factorisation.
pub(crate) aij: Vec<f64>,
/// Diagonal values, size `n`. Overwritten with L[j,j] during factorisation.
pub(crate) aii: Vec<f64>,
/// RHS on input, solution (permuted H) on output.
pub(crate) f: Vec<f64>,
// Working arrays for the factorisation kernel.
pub(super) link_chain: Vec<usize>,
pub(super) first_ptr: Vec<usize>,
pub(super) temp: Vec<f64>,
pub(super) last_timings: SparsePhaseTimings,
}
impl SparseSolver {
/// Constructs the solver for the given junction adjacency graph.
///
/// `adj[i]` = set of 0-based junction indices adjacent to junction i.
pub fn new(n: usize, adj: &[BTreeSet<usize>]) -> Self {
if n == 0 {
return SparseSolver {
n: 0,
n_coeff: 0,
row: vec![],
xlnz: vec![0],
nzsub: vec![],
aij: vec![],
aii: vec![],
f: vec![],
link_chain: vec![],
first_ptr: vec![],
temp: vec![],
last_timings: SparsePhaseTimings::default(),
};
}
let (order, row) = match genmmd_order(n, adj) {
Ok(result) => result,
Err(_) => greedy_mdo(n, adj),
};
let mut filled_adj = adj.to_vec();
symbolic_fill(n, &order, &mut filled_adj);
let (xlnz, nzsub, n_coeff) = build_csc(n, &order, &row, &filled_adj);
SparseSolver {
n,
n_coeff,
row,
xlnz,
nzsub,
aij: vec![0.0; n_coeff],
aii: vec![0.0; n],
f: vec![0.0; n],
link_chain: vec![n; n],
first_ptr: vec![0; n],
temp: vec![0.0; n],
last_timings: SparsePhaseTimings::default(),
}
}
/// Resets all matrix values to zero before assembly.
pub fn clear(&mut self) {
self.aij.fill(0.0);
self.aii.fill(0.0);
self.f.fill(0.0);
}
/// Builds a lookup map from `(permuted_col, permuted_row)` → position in `aij`.
///
/// Used at construction time to compute `link_aij_pos` in `SolverContext`.
#[allow(dead_code)]
pub fn pos_map(&self) -> HashMap<(usize, usize), usize> {
let mut map = HashMap::with_capacity(self.n_coeff);
for col in 0..self.n {
for pos in self.xlnz[col]..self.xlnz[col + 1] {
let r = self.nzsub[pos];
map.insert((col, r), pos);
}
}
map
}
/// Numerically factorises A = LLᵀ and solves Lz = F, Lᵀx = z in place (§3.6).
///
/// On entry: `aii` = diagonal of A, `aij` = off-diagonal entries (assembled as
/// negative values, matching EPANET's sign convention), `f` = RHS.
/// On success: `f` holds the solution (permuted head vector).
/// Returns `true` on success, `false` if the matrix is ill-conditioned.
pub fn factorize_solve(&mut self) -> bool {
if solve_timing_enabled() {
return self.factorize_solve_timed();
}
self.factorize_solve_fast()
}
fn factorize_solve_fast(&mut self) -> bool {
self.factorize_solve_scalar()
}
/// Scalar left-looking Cholesky factorisation + solve (§3.6 EPANET kernel).
fn factorize_solve_scalar(&mut self) -> bool {
let n = self.n;
self.link_chain.fill(n);
self.first_ptr.fill(0);
// `temp` is zeroed entry-by-entry during the column scatter loop below
// (each written entry is reset to 0.0 before moving to the next column).
// By the time this function returns, `temp` is all-zeros again, so the
// upfront fill is redundant on every call after the first.
debug_assert!(
self.temp.iter().all(|&v| v == 0.0),
"factorize_solve_scalar: temp buffer not clean on entry"
);
let xlnz = &self.xlnz;
let nzsub = &self.nzsub;
let aij = &mut self.aij;
let aii = &mut self.aii;
let f = &mut self.f;
let temp = &mut self.temp;
let link_chain = &mut self.link_chain;
let first_ptr = &mut self.first_ptr;
debug_assert!(self.nzsub.iter().all(|&s| s < n));
debug_assert!(*xlnz.last().unwrap_or(&0) <= aij.len());
for j in 0..n {
let mut diagj = 0.0f64;
let mut k = link_chain[j];
while k != n {
let newk = link_chain[k];
let kfirst = first_ptr[k];
let ljk = aij[kfirst];
diagj += ljk * ljk;
let istrt = kfirst + 1;
let istop = xlnz[k + 1];
if istrt < istop {
first_ptr[k] = istrt;
let isub = nzsub[istrt];
link_chain[k] = link_chain[isub];
link_chain[isub] = k;
// SAFETY: `istrt..istop` are within `aij` and `nzsub` (both
// length `n_coeff`, bounded by `xlnz[k+1]` ≤ n_coeff).
// `nzsub` entries are all < n, so `temp.get_unchecked_mut(isub)`
// is in bounds. The debug_assert! above validates both invariants.
unsafe {
let mut row_ptr = nzsub.as_ptr().add(istrt);
let mut val_ptr = aij.as_ptr().add(istrt);
let row_end = nzsub.as_ptr().add(istop);
while row_ptr != row_end {
let isub = *row_ptr;
*temp.get_unchecked_mut(isub) += *val_ptr * ljk;
row_ptr = row_ptr.add(1);
val_ptr = val_ptr.add(1);
}
}
}
k = newk;
}
diagj = aii[j] - diagj;
if diagj <= 0.0 {
return false;
}
diagj = diagj.sqrt();
aii[j] = diagj;
let istrt = xlnz[j];
let istop = xlnz[j + 1];
if istrt < istop {
first_ptr[j] = istrt;
let isub = nzsub[istrt];
link_chain[j] = link_chain[isub];
link_chain[isub] = j;
// SAFETY: same invariants as the first unsafe block above:
// `istrt..istop` ⊆ `0..n_coeff`, all `nzsub` entries < n,
// so all `temp` / `aij` / `f` accesses stay in bounds.
unsafe {
let mut row_ptr = nzsub.as_ptr().add(istrt);
let mut val_ptr = aij.as_mut_ptr().add(istrt);
let row_end = nzsub.as_ptr().add(istop);
while row_ptr != row_end {
let isub = *row_ptr;
let bj = (*val_ptr - *temp.get_unchecked(isub)) / diagj;
*val_ptr = bj;
*temp.get_unchecked_mut(isub) = 0.0;
row_ptr = row_ptr.add(1);
val_ptr = val_ptr.add(1);
}
}
}
}
for j in 0..n {
let bj = f[j] / aii[j];
f[j] = bj;
// SAFETY: `i` ranges over `xlnz[j]..xlnz[j+1]` ⊆ `0..n_coeff`;
// `nzsub[i] < n` (debug_assert! above) so `f.get_unchecked_mut` is valid.
unsafe {
let mut row_ptr = nzsub.as_ptr().add(xlnz[j]);
let mut val_ptr = aij.as_ptr().add(xlnz[j]);
let row_end = nzsub.as_ptr().add(xlnz[j + 1]);
while row_ptr != row_end {
let isub = *row_ptr;
*f.get_unchecked_mut(isub) -= *val_ptr * bj;
row_ptr = row_ptr.add(1);
val_ptr = val_ptr.add(1);
}
}
}
for j in (0..n).rev() {
let mut bj = f[j];
// SAFETY: same range invariants as the forward solve above;
// `f.get_unchecked(isub)` is valid because all `nzsub` entries < n.
unsafe {
let mut row_ptr = nzsub.as_ptr().add(xlnz[j]);
let mut val_ptr = aij.as_ptr().add(xlnz[j]);
let row_end = nzsub.as_ptr().add(xlnz[j + 1]);
while row_ptr != row_end {
let isub = *row_ptr;
bj -= *val_ptr * *f.get_unchecked(isub);
row_ptr = row_ptr.add(1);
val_ptr = val_ptr.add(1);
}
}
f[j] = bj / aii[j];
}
true
}
fn factorize_solve_timed(&mut self) -> bool {
let mut timings = SparsePhaseTimings::default();
let n = self.n;
let phase_started = Instant::now();
self.link_chain.fill(n);
self.first_ptr.fill(0);
self.temp.fill(0.0);
timings.reset += phase_started.elapsed();
debug_assert!(self.nzsub.iter().all(|&s| s < n));
debug_assert!(*self.xlnz.last().unwrap_or(&0) <= self.aij.len());
let phase_started = Instant::now();
for j in 0..n {
let mut diagj = 0.0f64;
let mut k = self.link_chain[j];
while k != n {
let newk = self.link_chain[k];
let kfirst = self.first_ptr[k];
let ljk = self.aij[kfirst];
diagj += ljk * ljk;
let istrt = kfirst + 1;
let istop = self.xlnz[k + 1];
if istrt < istop {
self.first_ptr[k] = istrt;
let isub = self.nzsub[istrt];
self.link_chain[k] = self.link_chain[isub];
self.link_chain[isub] = k;
// SAFETY: same as the factorisation inner loop above:
// all indices bounded by CSC structure invariants, nzsub entries < n. // SAFETY: `istrt..istop` ⊆ `0..n_coeff`; `nzsub` entries < n,
// so `temp.get_unchecked_mut(isub)` / `aij.get_unchecked(i)` are valid.
unsafe {
for i in istrt..istop {
let isub = *self.nzsub.get_unchecked(i);
*self.temp.get_unchecked_mut(isub) += *self.aij.get_unchecked(i) * ljk;
}
}
}
k = newk;
}
diagj = self.aii[j] - diagj;
if diagj <= 0.0 {
return false;
}
diagj = diagj.sqrt();
self.aii[j] = diagj;
let istrt = self.xlnz[j];
let istop = self.xlnz[j + 1];
if istrt < istop {
self.first_ptr[j] = istrt;
let isub = self.nzsub[istrt];
self.link_chain[j] = self.link_chain[isub];
self.link_chain[isub] = j;
unsafe {
for i in istrt..istop {
let isub = *self.nzsub.get_unchecked(i);
let bj =
(*self.aij.get_unchecked(i) - *self.temp.get_unchecked(isub)) / diagj;
*self.aij.get_unchecked_mut(i) = bj;
*self.temp.get_unchecked_mut(isub) = 0.0;
}
}
}
}
timings.factor += phase_started.elapsed();
let phase_started = Instant::now();
for j in 0..n {
let bj = self.f[j] / self.aii[j];
self.f[j] = bj;
// SAFETY: `i` ∈ `xlnz[j]..xlnz[j+1]` ⊆ `0..n_coeff`;
// `nzsub[i] < n` (debug_assert! above) so f/nzsub/aij accesses are valid.
unsafe {
for i in self.xlnz[j]..self.xlnz[j + 1] {
let isub = *self.nzsub.get_unchecked(i);
*self.f.get_unchecked_mut(isub) -= *self.aij.get_unchecked(i) * bj;
}
}
}
timings.forward += phase_started.elapsed();
let phase_started = Instant::now();
for j in (0..n).rev() {
let mut bj = self.f[j];
// SAFETY: same invariants as forward solve above.
unsafe {
for i in self.xlnz[j]..self.xlnz[j + 1] {
let isub = *self.nzsub.get_unchecked(i);
bj -= *self.aij.get_unchecked(i) * *self.f.get_unchecked(isub);
}
}
self.f[j] = bj / self.aii[j];
}
self.last_timings = timings;
self.last_timings.backward += phase_started.elapsed();
true
}
}
#[cfg(test)]
mod tests {
use super::*;
#[test]
fn sparse_solver_new_zero_system_has_empty_storage() {
let solver = SparseSolver::new(0, &[]);
assert_eq!(solver.n, 0);
assert_eq!(solver.n_coeff, 0);
assert_eq!(solver.xlnz, vec![0]);
assert!(solver.aij.is_empty());
assert!(solver.aii.is_empty());
}
#[test]
fn sparse_solver_clear_resets_matrix_and_rhs() {
let adj = vec![BTreeSet::from([1usize]), BTreeSet::from([0usize])];
let mut solver = SparseSolver::new(2, &adj);
solver.aij.fill(3.0);
solver.aii.fill(4.0);
solver.f.fill(5.0);
solver.clear();
assert!(solver.aij.iter().all(|v| *v == 0.0));
assert!(solver.aii.iter().all(|v| *v == 0.0));
assert!(solver.f.iter().all(|v| *v == 0.0));
}
#[test]
fn sparse_solver_pos_map_contains_lower_triangle_edge() {
let adj = vec![BTreeSet::from([1usize]), BTreeSet::from([0usize])];
let solver = SparseSolver::new(2, &adj);
let map = solver.pos_map();
assert_eq!(map.len(), 1);
assert!(map.contains_key(&(0, 1)) || map.contains_key(&(1, 0)));
}
}