ATLAS Offline Software
Loading...
Searching...
No Matches
ALFA_SvdCalc.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5/*
6 * svdcomp - SVD decomposition routine.
7 * Takes an mxn matrix a and decomposes it into udv, where u,v are
8 * left and right orthogonal transformation matrices, and d is a
9 * diagonal matrix of singular values.
10 *
11 * This routine is adapted from svdecomp.c in XLISP-STAT 2.1 which is
12 * code from Numerical Recipes adapted by Luke Tierney and David Betz.
13 *
14 * Input to dsvd is as follows:
15 * a = mxn matrix to be decomposed, gets overwritten with u
16 * m = row dimension of a
17 * n = column dimension of a
18 * w = returns the vector of singular values of a
19 * v = returns the right orthogonal transformation matrix
20*/
21
23#include <cmath>
24#include <cstdio>
25#include <cstdlib>
26#include <vector>
27
28static double PYTHAG(double a, double b)
29{
30 double at = fabs(a), bt = fabs(b), ct, result;
31
32 if (at > bt) { ct = bt / at; result = at * sqrt(1.0 + ct * ct); }
33 else if (bt > 0.0) { ct = at / bt; result = bt * sqrt(1.0 + ct * ct); }
34 else result = 0.0;
35 return(result);
36}
37
38
39int dsvd(double **a, int m, int n, double *w, double **v)
40{
41 int flag, i, its, j, jj, k, l=0, nm=0;
42 double c, f, h, s, x, y, z;
43 double anorm = 0.0, g = 0.0, scale = 0.0;
44
45 if (m < n)
46 {
47 fprintf(stderr, "#rows must be > #cols \n");
48 return(0);
49 }
50
51 std::vector<double> rv1 (n);
52
53/* Householder reduction to bidiagonal form */
54 for (i = 0; i < n; i++)
55 {
56 /* left-hand reduction */
57 l = i + 1;
58 rv1[i] = scale * g;
59 g = s = scale = 0.0;
60 if (i < m)
61 {
62 for (k = i; k < m; k++)
63 scale += fabs((double)a[k][i]);
64 if (scale)
65 {
66 for (k = i; k < m; k++)
67 {
68 a[k][i] = (double)((double)a[k][i]/scale);
69 s += ((double)a[k][i] * (double)a[k][i]);
70 }
71 f = (double)a[i][i];
72 g = -SIGN(sqrt(s), f);
73 h = f * g - s;
74 a[i][i] = (double)(f - g);
75 if (i != n - 1)
76 {
77 for (j = l; j < n; j++)
78 {
79 for (s = 0.0, k = i; k < m; k++)
80 s += ((double)a[k][i] * (double)a[k][j]);
81 f = s / h;
82 for (k = i; k < m; k++)
83 a[k][j] += (double)(f * (double)a[k][i]);
84 }
85 }
86 for (k = i; k < m; k++)
87 a[k][i] = (double)((double)a[k][i]*scale);
88 }
89 }
90 w[i] = (double)(scale * g);
91
92 /* right-hand reduction */
93 g = s = scale = 0.0;
94 if (i < m && i != n - 1)
95 {
96 for (k = l; k < n; k++)
97 scale += fabs((double)a[i][k]);
98 if (scale)
99 {
100 for (k = l; k < n; k++)
101 {
102 a[i][k] = (double)((double)a[i][k]/scale);
103 s += ((double)a[i][k] * (double)a[i][k]);
104 }
105 f = (double)a[i][l];
106 g = -SIGN(sqrt(s), f);
107 h = f * g - s;
108 a[i][l] = (double)(f - g);
109 for (k = l; k < n; k++)
110 rv1[k] = (double)a[i][k] / h;
111 if (i != m - 1)
112 {
113 for (j = l; j < m; j++)
114 {
115 for (s = 0.0, k = l; k < n; k++)
116 s += ((double)a[j][k] * (double)a[i][k]);
117 for (k = l; k < n; k++)
118 a[j][k] += (double)(s * rv1[k]);
119 }
120 }
121 for (k = l; k < n; k++)
122 a[i][k] = (double)((double)a[i][k]*scale);
123 }
124 }
125 anorm = MAX(anorm, (fabs((double)w[i]) + fabs(rv1[i])));
126 }
127
128 /* accumulate the right-hand transformation */
129 for (i = n - 1; i >= 0; i--)
130 {
131 if (i < n - 1)
132 {
133 if (g)
134 {
135 for (j = l; j < n; j++)
136 v[j][i] = (double)(((double)a[i][j] / (double)a[i][l]) / g);
137 /* double division to avoid underflow */
138 for (j = l; j < n; j++)
139 {
140 for (s = 0.0, k = l; k < n; k++)
141 s += ((double)a[i][k] * (double)v[k][j]);
142 for (k = l; k < n; k++)
143 v[k][j] += (double)(s * (double)v[k][i]);
144 }
145 }
146 for (j = l; j < n; j++)
147 v[i][j] = v[j][i] = 0.0;
148 }
149 v[i][i] = 1.0;
150 g = rv1[i];
151 l = i;
152 }
153
154 /* accumulate the left-hand transformation */
155 for (i = n - 1; i >= 0; i--)
156 {
157 l = i + 1;
158 g = (double)w[i];
159 if (i < n - 1)
160 for (j = l; j < n; j++)
161 a[i][j] = 0.0;
162 if (g)
163 {
164 g = 1.0 / g;
165 if (i != n - 1)
166 {
167 for (j = l; j < n; j++)
168 {
169 for (s = 0.0, k = l; k < m; k++)
170 s += ((double)a[k][i] * (double)a[k][j]);
171 f = (s / (double)a[i][i]) * g;
172 for (k = i; k < m; k++)
173 a[k][j] += (double)(f * (double)a[k][i]);
174 }
175 }
176 for (j = i; j < m; j++)
177 a[j][i] = (double)((double)a[j][i]*g);
178 }
179 else
180 {
181 for (j = i; j < m; j++)
182 a[j][i] = 0.0;
183 }
184 ++a[i][i];
185 }
186
187 /* diagonalize the bidiagonal form */
188 for (k = n - 1; k >= 0; k--)
189 { /* loop over singular values */
190 for (its = 0; its < 30; its++)
191 { /* loop over allowed iterations */
192 flag = 1;
193 for (l = k; l >= 0; l--)
194 { /* test for splitting */
195 nm = l - 1;
196 if (fabs(rv1[l]) + anorm == anorm)
197 {
198 flag = 0;
199 break;
200 }
201 if (fabs((double)w[nm]) + anorm == anorm)
202 break;
203 }
204 if (flag)
205 {
206 c = 0.0;
207 s = 1.0;
208 for (i = l; i <= k; i++)
209 {
210 f = s * rv1[i];
211 if (fabs(f) + anorm != anorm)
212 {
213 g = (double)w[i];
214 h = PYTHAG(f, g);
215 if (h ==0)[[unlikely]]{
216 fprintf(stderr, "Numerator h is zero! \n");
217 }
218 w[i] = (double)h;
219 h = 1.0 / h;
220 c = g * h;
221 s = (- f * h);
222 for (j = 0; j < m; j++)
223 {
224 y = (double)a[j][nm];
225 z = (double)a[j][i];
226 a[j][nm] = (double)(y * c + z * s);
227 a[j][i] = (double)(z * c - y * s);
228 }
229 }
230 }
231 }
232 z = (double)w[k];
233 if (l == k)
234 { /* convergence */
235 if (z < 0.0)
236 { /* make singular value nonnegative */
237 w[k] = (double)(-z);
238 for (j = 0; j < n; j++)
239 v[j][k] = (-v[j][k]);
240 }
241 break;
242 }
243 if (its >= 30) {
244 fprintf(stderr, "No convergence after 30,000! iterations \n");
245 return(0);
246 }
247
248 /* shift from bottom 2 x 2 minor */
249 x = (double)w[l];
250 nm = k - 1;
251 y = (double)w[nm];
252 g = rv1[nm];
253 h = rv1[k];
254 f = ((y - z) * (y + z) + (g - h) * (g + h)) / (2.0 * h * y);
255 g = PYTHAG(f, 1.0);
256 f = ((x - z) * (x + z) + h * ((y / (f + SIGN(g, f))) - h)) / x;
257
258 /* next QR transformation */
259 c = s = 1.0;
260 for (j = l; j <= nm; j++)
261 {
262 i = j + 1;
263 g = rv1[i];
264 y = (double)w[i];
265 h = s * g;
266 g = c * g;
267 z = PYTHAG(f, h);
268 if (z == 0)[[unlikely]]{
269 fprintf(stderr, "Numerator z is zero! \n");
270 return(0);
271 }
272 rv1[j] = z;
273 c = f / z;
274 s = h / z;
275 f = x * c + g * s;
276 g = g * c - x * s;
277 h = y * s;
278 y = y * c;
279 for (jj = 0; jj < n; jj++)
280 {
281 x = (double)v[jj][j];
282 z = (double)v[jj][i];
283 v[jj][j] = (double)(x * c + z * s);
284 v[jj][i] = (double)(z * c - x * s);
285 }
286 z = PYTHAG(f, h);
287 w[j] = (double)z;
288 if (z)
289 {
290 z = 1.0 / z;
291 c = f * z;
292 s = h * z;
293 }
294 f = (c * g) + (s * y);
295 x = (c * y) - (s * g);
296 for (jj = 0; jj < m; jj++)
297 {
298 y = (double)a[jj][j];
299 z = (double)a[jj][i];
300 a[jj][j] = (double)(y * c + z * s);
301 a[jj][i] = (double)(z * c - y * s);
302 }
303 }
304 rv1[l] = 0.0;
305 rv1[k] = f;
306 w[k] = (double)x;
307 }
308 }
309 return(1);
310}
int dsvd(double **a, int m, int n, double *w, double **v)
static double PYTHAG(double a, double b)
#define MAX(x, y)
#define SIGN(a, b)
static Double_t a
#define y
#define x
#define z
Header file for AthHistogramAlgorithm.
#define unlikely(x)
Definition dictionary.h:41