source: git/src/matrix.c @ e3f7ea7

stereo-2025
Last change on this file since e3f7ea7 was 4e7fb5e, checked in by Olly Betts <olly@…>, 7 months ago

Improve comment about Choleski

  • Property mode set to 100644
File size: 12.7 KB
Line 
1/* matrix.c
2 * Matrix building and solving routines
3 * Copyright (C) 1993-2003,2010,2013,2024 Olly Betts
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
18 */
19
20/*#define SOR 1*/
21
22#if 0
23# define DEBUG_INVALID 1
24#endif
25
26#include <config.h>
27
28#include "debug.h"
29#include "cavern.h"
30#include "filename.h"
31#include "message.h"
32#include "netbits.h"
33#include "matrix.h"
34#include "out.h"
35
36#undef PRINT_MATRICES
37#define PRINT_MATRICES 0
38
39#undef DEBUG_MATRIX_BUILD
40#define DEBUG_MATRIX_BUILD 0
41
42#undef DEBUG_MATRIX
43#define DEBUG_MATRIX 0
44
45#if PRINT_MATRICES
46static void print_matrix(real *M, real *B, long n);
47#endif
48
49static void choleski(real *M, real *B, long n);
50
51#ifdef SOR
52static void sor(real *M, real *B, long n);
53#endif
54
55/* for M(row, col) col must be <= row, so Y <= X */
56# define M(X, Y) ((real *)M)[((((OSSIZE_T)(X)) * ((X) + 1)) >> 1) + (Y)]
57              /* +(Y>X?0*printf("row<col (line %d)\n",__LINE__):0) */
58/*#define M_(X, Y) ((real *)M)[((((OSSIZE_T)(Y)) * ((Y) + 1)) >> 1) + (X)]*/
59
60static void set_row(node *stn, int row_number) {
61    // We store the matrix row/column index in stn->colour for quick and easy
62    // lookup when copying out the solved station coordinates.
63    stn->colour = row_number;
64    for (int d = 0; d < 3; d++) {
65        linkfor *leg = stn->leg[d];
66        if (!leg) break;
67        node *to = leg->l.to;
68        if (to->colour < 0 && stn->name->pos == to->name->pos) {
69            set_row(to, row_number);
70        }
71    }
72}
73
74#ifdef NO_COVARIANCES
75# define FACTOR 1
76#else
77# define FACTOR 3
78#endif
79
80/* Find positions for a subset of the reduced network by solving a matrix
81 * equation.
82 *
83 * list is a non-empty linked list of unfixed stations to solve for.
84 *
85 * As a pre-condition, all stations in list must have a negative value for
86 * stn->colour.  This can be ensured by the caller (which avoids having to
87 * make an extra pass over the list just to set the colours suitably).
88 */
89extern void
90solve_matrix(node *list)
91{
92   // Assign a matrix row/column index to each group of stations with the same
93   // pos.
94   //
95   // We also set listend to the last station in the list while doing so, which
96   // we use after solving to splice list into fixedlist.
97   node *listend = NULL;
98   long n = 0;
99   for (node *stn = list; stn; stn = stn->next) {
100      listend = stn;
101      if (stn->colour < 0) {
102          set_row(stn, n++);
103      }
104   }
105   SVX_ASSERT(n > 0);
106
107   // Array to map from row/column index to pos.  We fill this in as we build
108   // the matrix, and use it to know where to copy the solved station
109   // coordinates to.
110   pos **stn_tab = osmalloc((OSSIZE_T)(n * ossizeof(pos*)));
111
112   /* (OSSIZE_T) cast may be needed if n >= 181 */
113   real *M = osmalloc((OSSIZE_T)((((OSSIZE_T)n * FACTOR * (n * FACTOR + 1)) >> 1)) * ossizeof(real));
114   real *B = osmalloc((OSSIZE_T)(n * FACTOR * ossizeof(real)));
115
116   if (!fQuiet) {
117      if (n == 1)
118         out_current_action(msg(/*Solving one equation*/78));
119      else
120         out_current_action1(msg(/*Solving %d simultaneous equations*/75), n);
121   }
122
123#ifdef NO_COVARIANCES
124   int dim = 2;
125#else
126   int dim = 0; /* Collapse loop to a single iteration. */
127#endif
128   for ( ; dim >= 0; dim--) {
129      /* Initialise M and B to zero - zeroing "linearly" will minimise
130       * paging when the matrix is large */
131      {
132         int end = n * FACTOR;
133         for (int row = 0; row < end; row++) B[row] = (real)0.0;
134         end = ((OSSIZE_T)n * FACTOR * (n * FACTOR + 1)) >> 1;
135         for (int row = 0; row < end; row++) M[row] = (real)0.0;
136      }
137
138      /* Construct matrix by going through the stn list.
139       *
140       * All legs between two fixed stations can be ignored here.
141       *
142       * Other legs we want to add exactly once to M.  To achieve this we
143       * want to:
144       *
145       * - add forward legs between two unfixed stations,
146       *
147       * - add legs from unfixed stations to fixed stations (we do them from
148       *   the unfixed end so we don't need to detect when we're at a fixed
149       *   point cut line and determine which side we're currently dealing
150       *   with).
151       *
152       * To implement this, we only look at legs from unfixed stations and add
153       * a leg if to a fixed station, or to an unfixed station and it's a
154       * forward leg.
155       */
156      for (node *stn = list; stn; stn = stn->next) {
157         if (dim == 0) {
158             stn_tab[stn->colour] = stn->name->pos;
159         }
160
161#ifdef NO_COVARIANCES
162         real e;
163#else
164         svar e;
165         delta a;
166#endif
167#if DEBUG_MATRIX_BUILD
168         print_prefix(stn->name);
169         printf(" used: %d colour %ld\n",
170                (!!stn->leg[2]) << 2 | (!!stn -> leg[1]) << 1 | (!!stn->leg[0]),
171                stn->colour);
172
173         for (int dirn = 0; dirn <= 2 && stn->leg[dirn]; dirn++) {
174            printf("Leg %d, vx=%f, reverse=%d, to ", dirn,
175                   stn->leg[dirn]->v[0], stn->leg[dirn]->l.reverse);
176            print_prefix(stn->leg[dirn]->l.to->name);
177            putnl();
178         }
179         putnl();
180#endif /* DEBUG_MATRIX_BUILD */
181
182         int f = stn->colour;
183         SVX_ASSERT(f >= 0);
184         {
185            for (int dirn = 0; dirn <= 2 && stn->leg[dirn]; dirn++) {
186               linkfor *leg = stn->leg[dirn];
187               node *to = leg->l.to;
188               if (fixed(to)) {
189                  bool fRev = !data_here(leg);
190                  if (fRev) leg = reverse_leg(leg);
191                  /* Ignore equated nodes */
192#ifdef NO_COVARIANCES
193                  e = leg->v[dim];
194                  if (e != (real)0.0) {
195                     e = ((real)1.0) / e;
196                     M(f,f) += e;
197                     B[f] += e * POS(to, dim);
198                     if (fRev) {
199                        B[f] += leg->d[dim];
200                     } else {
201                        B[f] -= leg->d[dim];
202                     }
203                  }
204#else
205                  if (invert_svar(&e, &leg->v)) {
206                     if (fRev) {
207                        adddd(&a, &POSD(to), &leg->d);
208                     } else {
209                        subdd(&a, &POSD(to), &leg->d);
210                     }
211                     delta b;
212                     mulsd(&b, &e, &a);
213                     for (int i = 0; i < 3; i++) {
214                        M(f * FACTOR + i, f * FACTOR + i) += e[i];
215                        B[f * FACTOR + i] += b[i];
216                     }
217                     M(f * FACTOR + 1, f * FACTOR) += e[3];
218                     M(f * FACTOR + 2, f * FACTOR) += e[4];
219                     M(f * FACTOR + 2, f * FACTOR + 1) += e[5];
220                  }
221#endif
222               } else if (data_here(leg) &&
223                          (leg->l.reverse & FLAG_ARTICULATION) == 0) {
224                  /* forward leg, unfixed -> unfixed */
225                  int t = to->colour;
226                  SVX_ASSERT(t >= 0);
227#if DEBUG_MATRIX
228# ifdef NO_COVARIANCES
229                  printf("Leg %d to %d, var %f, delta %f\n", f, t, e,
230                         leg->d[dim]);
231# else
232                  printf("Leg %d to %d, var (%f, %f, %f; %f, %f, %f), "
233                         "delta %f\n", f, t, e[0], e[1], e[2], e[3], e[4], e[5],
234                         leg->d[dim]);
235# endif
236#endif
237                  /* Ignore equated nodes & lollipops */
238#ifdef NO_COVARIANCES
239                  e = leg->v[dim];
240                  if (t != f && e != (real)0.0) {
241                     e = ((real)1.0) / e;
242                     M(f,f) += e;
243                     M(t,t) += e;
244                     if (f < t) M(t,f) -= e; else M(f,t) -= e;
245                     real a = e * leg->d[dim];
246                     B[f] -= a;
247                     B[t] += a;
248                  }
249#else
250                  if (t != f && invert_svar(&e, &leg->v)) {
251                     mulsd(&a, &e, &leg->d);
252                     for (int i = 0; i < 3; i++) {
253                        M(f * FACTOR + i, f * FACTOR + i) += e[i];
254                        M(t * FACTOR + i, t * FACTOR + i) += e[i];
255                        if (f < t)
256                           M(t * FACTOR + i, f * FACTOR + i) -= e[i];
257                        else
258                           M(f * FACTOR + i, t * FACTOR + i) -= e[i];
259                        B[f * FACTOR + i] -= a[i];
260                        B[t * FACTOR + i] += a[i];
261                     }
262                     M(f * FACTOR + 1, f * FACTOR) += e[3];
263                     M(t * FACTOR + 1, t * FACTOR) += e[3];
264                     M(f * FACTOR + 2, f * FACTOR) += e[4];
265                     M(t * FACTOR + 2, t * FACTOR) += e[4];
266                     M(f * FACTOR + 2, f * FACTOR + 1) += e[5];
267                     M(t * FACTOR + 2, t * FACTOR + 1) += e[5];
268                     if (f < t) {
269                        M(t * FACTOR + 1, f * FACTOR) -= e[3];
270                        M(t * FACTOR, f * FACTOR + 1) -= e[3];
271                        M(t * FACTOR + 2, f * FACTOR) -= e[4];
272                        M(t * FACTOR, f * FACTOR + 2) -= e[4];
273                        M(t * FACTOR + 2, f * FACTOR + 1) -= e[5];
274                        M(t * FACTOR + 1, f * FACTOR + 2) -= e[5];
275                     } else {
276                        M(f * FACTOR + 1, t * FACTOR) -= e[3];
277                        M(f * FACTOR, t * FACTOR + 1) -= e[3];
278                        M(f * FACTOR + 2, t * FACTOR) -= e[4];
279                        M(f * FACTOR, t * FACTOR + 2) -= e[4];
280                        M(f * FACTOR + 2, t * FACTOR + 1) -= e[5];
281                        M(f * FACTOR + 1, t * FACTOR + 2) -= e[5];
282                     }
283                  }
284#endif
285               }
286            }
287         }
288      }
289
290#if PRINT_MATRICES
291      print_matrix(M, B, n * FACTOR); /* 'ave a look! */
292#endif
293
294#ifdef SOR
295      /* defined in network.c, may be altered by -z<letters> on command line */
296      if (optimize & BITA('i'))
297         sor(M, B, n * FACTOR);
298      else
299#endif
300         choleski(M, B, n * FACTOR);
301
302      {
303         for (int m = (int)(n - 1); m >= 0; m--) {
304#ifdef NO_COVARIANCES
305            stn_tab[m]->p[dim] = B[m];
306            if (dim == 0) {
307               SVX_ASSERT2(pos_fixed(stn_tab[m]),
308                       "setting station coordinates didn't mark pos as fixed");
309            }
310#else
311            for (int i = 0; i < 3; i++) {
312               stn_tab[m]->p[i] = B[m * FACTOR + i];
313            }
314            SVX_ASSERT2(pos_fixed(stn_tab[m]),
315                    "setting station coordinates didn't mark pos as fixed");
316#endif
317         }
318      }
319   }
320
321   // Put the solved stations back on fixedlist.
322   listend->next = fixedlist;
323   if (fixedlist) fixedlist->prev = listend;
324   fixedlist = list;
325
326   osfree(B);
327   osfree(M);
328   osfree(stn_tab);
329
330#if DEBUG_MATRIX
331   for (node *stn = list; stn; stn = stn->next) {
332      printf("(%8.2f, %8.2f, %8.2f ) ", POS(stn, 0), POS(stn, 1), POS(stn, 2));
333      print_prefix(stn->name);
334      putnl();
335   }
336#endif
337}
338
339/* Solve MX=B for X by first factoring M into LDL'.  This is a modified form
340 * of Choleski factorisation - the original Choleski factorisation is LL',
341 * but this modified version has the advantage of avoiding O(n) square root
342 * calculations.
343 */
344/* Note M must be symmetric positive definite */
345/* routine is entitled to scribble on M and B if it wishes */
346static void
347choleski(real *M, real *B, long n)
348{
349   for (int j = 1; j < n; j++) {
350      real V;
351      for (int i = 0; i < j; i++) {
352         V = (real)0.0;
353         for (int k = 0; k < i; k++) V += M(i,k) * M(j,k) * M(k,k);
354         M(j,i) = (M(j,i) - V) / M(i,i);
355      }
356      V = (real)0.0;
357      for (int k = 0; k < j; k++) V += M(j,k) * M(j,k) * M(k,k);
358      M(j,j) -= V; /* may be best to add M() last for numerical reasons too */
359   }
360
361   /* Multiply x by L inverse */
362   for (int i = 0; i < n - 1; i++) {
363      for (int j = i + 1; j < n; j++) {
364         B[j] -= M(j,i) * B[i];
365      }
366   }
367
368   /* Multiply x by D inverse */
369   for (int i = 0; i < n; i++) {
370      B[i] /= M(i,i);
371   }
372
373   /* Multiply x by (L transpose) inverse */
374   for (int i = (int)(n - 1); i > 0; i--) {
375      for (int j = i - 1; j >= 0; j--) {
376         B[j] -= M(i,j) * B[i];
377      }
378   }
379
380   /* printf("\n%ld/%ld\n\n",flops,flopsTot); */
381}
382
383#ifdef SOR
384/* factor to use for SOR (must have 1 <= SOR_factor < 2) */
385#define SOR_factor 1.93 /* 1.95 */
386
387/* Solve MX=B for X by SOR of Gauss-Siedel */
388/* routine is entitled to scribble on M and B if it wishes */
389static void
390sor(real *M, real *B, long n)
391{
392   long it = 0;
393
394   real *X = osmalloc(n * ossizeof(real));
395
396   const real threshold = 0.00001;
397
398   printf("reciprocating diagonal\n"); /* TRANSLATE */
399
400   /* munge diagonal so we can multiply rather than divide */
401   for (int row = n - 1; row >= 0; row--) {
402      M(row,row) = 1 / M(row,row);
403      X[row] = 0;
404   }
405
406   printf("starting iteration\n"); /* TRANSLATE */
407
408   real t;
409   do {
410      /*printf("*");*/
411      it++;
412      t = 0.0;
413      for (int row = 0; row < n; row++) {
414         real x = B[row];
415         int col;
416         for (col = 0; col < row; col++) x -= M(row,col) * X[col];
417         for (col++; col < n; col++) x -= M(col,row) * X[col];
418         x *= M(row,row);
419         real sor_delta = (x - X[row]) * SOR_factor;
420         X[row] += sor_delta;
421         real t2 = fabs(sor_delta);
422         if (t2 > t) t = t2;
423      }
424      printf("% 6ld: %8.6f\n", it, t);
425   } while (t >= threshold && it < 100000);
426
427   if (t >= threshold) {
428      fprintf(stderr, "*not* converged after %ld iterations\n", it);
429      BUG("iteration stinks");
430   }
431
432   printf("%ld iterations\n", it); /* TRANSLATE */
433
434#if 0
435   putnl();
436   for (int row = n - 1; row >= 0; row--) {
437      t = 0.0;
438      for (int col = 0; col < row; col++) t += M(row, col) * X[col];
439      t += X[row] / M(row, row);
440      for (col = row + 1; col < n; col++)
441         t += M(col, row) * X[col];
442      printf("[ %f %f ]\n", t, B[row]);
443   }
444#endif
445
446   for (int row = n - 1; row >= 0; row--) B[row] = X[row];
447
448   osfree(X);
449   printf("\ndone\n"); /* TRANSLATE */
450}
451#endif
452
453#if PRINT_MATRICES
454static void
455print_matrix(real *M, real *B, long n)
456{
457   printf("Matrix, M and vector, B:\n");
458   for (long row = 0; row < n; row++) {
459      long col;
460      for (col = 0; col <= row; col++) printf("%6.2f\t", M(row, col));
461      for (; col <= n; col++) printf(" \t");
462      printf("\t%6.2f\n", B[row]);
463   }
464   putnl();
465   return;
466}
467#endif
Note: See TracBrowser for help on using the repository browser.