source: git/src/matrix.c @ 2d8d46d

stereo-2025
Last change on this file since 2d8d46d was 2d8d46d, checked in by Olly Betts <olly@…>, 9 months ago

Assign matrix rows more efficiently

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