NEURON
qrfactor.c
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 
3 /**************************************************************************
4 **
5 ** Copyright (C) 1993 David E. Steward & Zbigniew Leyk, all rights reserved.
6 **
7 ** Meschach Library
8 **
9 ** This Meschach Library is provided "as is" without any express
10 ** or implied warranty of any kind with respect to this software.
11 ** In particular the authors shall not be liable for any direct,
12 ** indirect, special, incidental or consequential damages arising
13 ** in any way from use of the software.
14 **
15 ** Everyone is granted permission to copy, modify and redistribute this
16 ** Meschach Library, provided:
17 ** 1. All copies contain this copyright notice.
18 ** 2. All modified copies shall carry a notice stating who
19 ** made the last modification and the date of such modification.
20 ** 3. No charge is made for this software or works derived from it.
21 ** This clause shall not be construed as constraining other software
22 ** distributed on the same medium as this software, nor is a
23 ** distribution fee considered a charge.
24 **
25 ***************************************************************************/
26 
27 
28 /*
29  This file contains the routines needed to perform QR factorisation
30  of matrices, as well as Householder transformations.
31  The internal "factored form" of a matrix A is not quite standard.
32  The diagonal of A is replaced by the diagonal of R -- not by the 1st non-zero
33  entries of the Householder vectors. The 1st non-zero entries are held in
34  the diag parameter of QRfactor(). The reason for this non-standard
35  representation is that it enables direct use of the Usolve() function
36  rather than requiring that a seperate function be written just for this case.
37  See, e.g., QRsolve() below for more details.
38 
39 */
40 
41 
42 static char rcsid[] = "qrfactor.c,v 1.1 1997/12/04 17:55:45 hines Exp";
43 
44 #include <stdio.h>
45 #include "matrix2.h"
46 #include <math.h>
47 
48 
49 
50 
51 
52 #define sign(x) ((x) > 0.0 ? 1 : ((x) < 0.0 ? -1 : 0 ))
53 
54 extern VEC *Usolve(); /* See matrix2.h */
55 
56 /* Note: The usual representation of a Householder transformation is taken
57  to be:
58  P = I - beta.u.uT
59  where beta = 2/(uT.u) and u is called the Householder vector
60  */
61 
62 /* QRfactor -- forms the QR factorisation of A -- factorisation stored in
63  compact form as described above ( not quite standard format ) */
64 /* MAT *QRfactor(A,diag,beta) */
66 MAT *A;
67 VEC *diag /* ,*beta */;
68 {
69  u_int k,limit;
70  Real beta;
71  static VEC *tmp1=VNULL;
72 
73  if ( ! A || ! diag )
74  error(E_NULL,"QRfactor");
75  limit = min(A->m,A->n);
76  if ( diag->dim < limit )
77  error(E_SIZES,"QRfactor");
78 
79  tmp1 = v_resize(tmp1,A->m);
80  MEM_STAT_REG(tmp1,TYPE_VEC);
81 
82  for ( k=0; k<limit; k++ )
83  {
84  /* get H/holder vector for the k-th column */
85  get_col(A,k,tmp1);
86  /* hhvec(tmp1,k,&beta->ve[k],tmp1,&A->me[k][k]); */
87  hhvec(tmp1,k,&beta,tmp1,&A->me[k][k]);
88  diag->ve[k] = tmp1->ve[k];
89 
90  /* apply H/holder vector to remaining columns */
91  /* hhtrcols(A,k,k+1,tmp1,beta->ve[k]); */
92  hhtrcols(A,k,k+1,tmp1,beta);
93  }
94 
95  return (A);
96 }
97 
98 /* QRCPfactor -- forms the QR factorisation of A with column pivoting
99  -- factorisation stored in compact form as described above
100  ( not quite standard format ) */
101 /* MAT *QRCPfactor(A,diag,beta,px) */
103 MAT *A;
104 VEC *diag /* , *beta */;
105 PERM *px;
106 {
107  u_int i, i_max, j, k, limit;
108  static VEC *gamma=VNULL, *tmp1=VNULL, *tmp2=VNULL;
109  Real beta, maxgamma, sum, tmp;
110 
111  if ( ! A || ! diag || ! px )
112  error(E_NULL,"QRCPfactor");
113  limit = min(A->m,A->n);
114  if ( diag->dim < limit || px->size != A->n )
115  error(E_SIZES,"QRCPfactor");
116 
117  tmp1 = v_resize(tmp1,A->m);
118  tmp2 = v_resize(tmp2,A->m);
119  gamma = v_resize(gamma,A->n);
120  MEM_STAT_REG(tmp1,TYPE_VEC);
121  MEM_STAT_REG(tmp2,TYPE_VEC);
122  MEM_STAT_REG(gamma,TYPE_VEC);
123 
124  /* initialise gamma and px */
125  for ( j=0; j<A->n; j++ )
126  {
127  px->pe[j] = j;
128  sum = 0.0;
129  for ( i=0; i<A->m; i++ )
130  sum += square(A->me[i][j]);
131  gamma->ve[j] = sum;
132  }
133 
134  for ( k=0; k<limit; k++ )
135  {
136  /* find "best" column to use */
137  i_max = k; maxgamma = gamma->ve[k];
138  for ( i=k+1; i<A->n; i++ )
139  /* Loop invariant:maxgamma=gamma[i_max]
140  >=gamma[l];l=k,...,i-1 */
141  if ( gamma->ve[i] > maxgamma )
142  { maxgamma = gamma->ve[i]; i_max = i; }
143 
144  /* swap columns if necessary */
145  if ( i_max != k )
146  {
147  /* swap gamma values */
148  tmp = gamma->ve[k];
149  gamma->ve[k] = gamma->ve[i_max];
150  gamma->ve[i_max] = tmp;
151 
152  /* update column permutation */
153  px_transp(px,k,i_max);
154 
155  /* swap columns of A */
156  for ( i=0; i<A->m; i++ )
157  {
158  tmp = A->me[i][k];
159  A->me[i][k] = A->me[i][i_max];
160  A->me[i][i_max] = tmp;
161  }
162  }
163 
164  /* get H/holder vector for the k-th column */
165  get_col(A,k,tmp1);
166  /* hhvec(tmp1,k,&beta->ve[k],tmp1,&A->me[k][k]); */
167  hhvec(tmp1,k,&beta,tmp1,&A->me[k][k]);
168  diag->ve[k] = tmp1->ve[k];
169 
170  /* apply H/holder vector to remaining columns */
171  /* hhtrcols(A,k,k+1,tmp1,beta->ve[k]); */
172  hhtrcols(A,k,k+1,tmp1,beta);
173 
174  /* update gamma values */
175  for ( j=k+1; j<A->n; j++ )
176  gamma->ve[j] -= square(A->me[k][j]);
177  }
178 
179  return (A);
180 }
181 
182 /* Qsolve -- solves Qx = b, Q is an orthogonal matrix stored in compact
183  form a la QRfactor() -- may be in-situ */
184 /* VEC *_Qsolve(QR,diag,beta,b,x,tmp) */
185 VEC *_Qsolve(QR,diag,b,x,tmp)
186 MAT *QR;
187 VEC *diag /* ,*beta */ , *b, *x, *tmp;
188 {
189  u_int dynamic;
190  int k, limit;
191  Real beta, r_ii, tmp_val;
192 
193  limit = min(QR->m,QR->n);
194  dynamic = FALSE;
195  if ( ! QR || ! diag || ! b )
196  error(E_NULL,"_Qsolve");
197  if ( diag->dim < limit || b->dim != QR->m )
198  error(E_SIZES,"_Qsolve");
199  x = v_resize(x,QR->m);
200  if ( tmp == VNULL )
201  dynamic = TRUE;
202  tmp = v_resize(tmp,QR->m);
203 
204  /* apply H/holder transforms in normal order */
205  x = v_copy(b,x);
206  for ( k = 0 ; k < limit ; k++ )
207  {
208  get_col(QR,k,tmp);
209  r_ii = fabs(tmp->ve[k]);
210  tmp->ve[k] = diag->ve[k];
211  tmp_val = (r_ii*fabs(diag->ve[k]));
212  beta = ( tmp_val == 0.0 ) ? 0.0 : 1.0/tmp_val;
213  /* hhtrvec(tmp,beta->ve[k],k,x,x); */
214  hhtrvec(tmp,beta,k,x,x);
215  }
216 
217  if ( dynamic )
218  V_FREE(tmp);
219 
220  return (x);
221 }
222 
223 /* makeQ -- constructs orthogonal matrix from Householder vectors stored in
224  compact QR form */
225 /* MAT *makeQ(QR,diag,beta,Qout) */
226 MAT *makeQ(QR,diag,Qout)
227 MAT *QR,*Qout;
228 VEC *diag /* , *beta */;
229 {
230  static VEC *tmp1=VNULL,*tmp2=VNULL;
231  u_int i, limit;
232  Real beta, r_ii, tmp_val;
233  int j;
234 
235  limit = min(QR->m,QR->n);
236  if ( ! QR || ! diag )
237  error(E_NULL,"makeQ");
238  if ( diag->dim < limit )
239  error(E_SIZES,"makeQ");
240  if ( Qout==(MAT *)NULL || Qout->m < QR->m || Qout->n < QR->m )
241  Qout = m_get(QR->m,QR->m);
242 
243  tmp1 = v_resize(tmp1,QR->m); /* contains basis vec & columns of Q */
244  tmp2 = v_resize(tmp2,QR->m); /* contains H/holder vectors */
245  MEM_STAT_REG(tmp1,TYPE_VEC);
246  MEM_STAT_REG(tmp2,TYPE_VEC);
247 
248  for ( i=0; i<QR->m ; i++ )
249  { /* get i-th column of Q */
250  /* set up tmp1 as i-th basis vector */
251  for ( j=0; j<QR->m ; j++ )
252  tmp1->ve[j] = 0.0;
253  tmp1->ve[i] = 1.0;
254 
255  /* apply H/h transforms in reverse order */
256  for ( j=limit-1; j>=0; j-- )
257  {
258  get_col(QR,j,tmp2);
259  r_ii = fabs(tmp2->ve[j]);
260  tmp2->ve[j] = diag->ve[j];
261  tmp_val = (r_ii*fabs(diag->ve[j]));
262  beta = ( tmp_val == 0.0 ) ? 0.0 : 1.0/tmp_val;
263  /* hhtrvec(tmp2,beta->ve[j],j,tmp1,tmp1); */
264  hhtrvec(tmp2,beta,j,tmp1,tmp1);
265  }
266 
267  /* insert into Q */
268  set_col(Qout,i,tmp1);
269  }
270 
271  return (Qout);
272 }
273 
274 /* makeR -- constructs upper triangular matrix from QR (compact form)
275  -- may be in-situ (all it does is zero the lower 1/2) */
276 MAT *makeR(QR,Rout)
277 MAT *QR,*Rout;
278 {
279  u_int i,j;
280 
281  if ( QR==(MAT *)NULL )
282  error(E_NULL,"makeR");
283  Rout = m_copy(QR,Rout);
284 
285  for ( i=1; i<QR->m; i++ )
286  for ( j=0; j<QR->n && j<i; j++ )
287  Rout->me[i][j] = 0.0;
288 
289  return (Rout);
290 }
291 
292 /* QRsolve -- solves the system Q.R.x=b where Q & R are stored in compact form
293  -- returns x, which is created if necessary */
294 /* VEC *QRsolve(QR,diag,beta,b,x) */
295 VEC *QRsolve(QR,diag,b,x)
296 MAT *QR;
297 VEC *diag /* , *beta */ , *b, *x;
298 {
299  int limit;
300  static VEC *tmp = VNULL;
301 
302  if ( ! QR || ! diag || ! b )
303  error(E_NULL,"QRsolve");
304  limit = min(QR->m,QR->n);
305  if ( diag->dim < limit || b->dim != QR->m )
306  error(E_SIZES,"QRsolve");
307  tmp = v_resize(tmp,limit);
308  MEM_STAT_REG(tmp,TYPE_VEC);
309 
310  x = v_resize(x,QR->n);
311  _Qsolve(QR,diag,b,x,tmp);
312  x = Usolve(QR,x,x,0.0);
313  v_resize(x,QR->n);
314 
315  return x;
316 }
317 
318 /* QRCPsolve -- solves A.x = b where A is factored by QRCPfactor()
319  -- assumes that A is in the compact factored form */
320 /* VEC *QRCPsolve(QR,diag,beta,pivot,b,x) */
321 VEC *QRCPsolve(QR,diag,pivot,b,x)
322 MAT *QR;
323 VEC *diag /* , *beta */;
324 PERM *pivot;
325 VEC *b, *x;
326 {
327  static VEC *tmp=VNULL;
328 
329  if ( ! QR || ! diag || ! pivot || ! b )
330  error(E_NULL,"QRCPsolve");
331  if ( (QR->m > diag->dim &&QR->n > diag->dim) || QR->n != pivot->size )
332  error(E_SIZES,"QRCPsolve");
333 
334  tmp = QRsolve(QR,diag /* , beta */ ,b,tmp);
335  MEM_STAT_REG(tmp,TYPE_VEC);
336  x = pxinv_vec(pivot,tmp,x);
337 
338  return x;
339 }
340 
341 /* Umlt -- compute out = upper_triang(U).x
342  -- may be in situ */
343 static VEC *Umlt(U,x,out)
344 MAT *U;
345 VEC *x, *out;
346 {
347  int i, limit;
348 
349  if ( U == MNULL || x == VNULL )
350  error(E_NULL,"Umlt");
351  limit = min(U->m,U->n);
352  if ( limit != x->dim )
353  error(E_SIZES,"Umlt");
354  if ( out == VNULL || out->dim < limit )
355  out = v_resize(out,limit);
356 
357  for ( i = 0; i < limit; i++ )
358  out->ve[i] = __ip__(&(x->ve[i]),&(U->me[i][i]),limit - i);
359  return out;
360 }
361 
362 /* UTmlt -- returns out = upper_triang(U)^T.x */
363 static VEC *UTmlt(U,x,out)
364 MAT *U;
365 VEC *x, *out;
366 {
367  Real sum;
368  int i, j, limit;
369 
370  if ( U == MNULL || x == VNULL )
371  error(E_NULL,"UTmlt");
372  limit = min(U->m,U->n);
373  if ( out == VNULL || out->dim < limit )
374  out = v_resize(out,limit);
375 
376  for ( i = limit-1; i >= 0; i-- )
377  {
378  sum = 0.0;
379  for ( j = 0; j <= i; j++ )
380  sum += U->me[j][i]*x->ve[j];
381  out->ve[i] = sum;
382  }
383  return out;
384 }
385 
386 /* QRTsolve -- solve A^T.sc = c where the QR factors of A are stored in
387  compact form
388  -- returns sc
389  -- original due to Mike Osborne modified Wed 09th Dec 1992 */
391 MAT *A;
392 VEC *diag, *c, *sc;
393 {
394  int i, j, k, n, p;
395  Real beta, r_ii, s, tmp_val;
396 
397  if ( ! A || ! diag || ! c )
398  error(E_NULL,"QRTsolve");
399  if ( diag->dim < min(A->m,A->n) )
400  error(E_SIZES,"QRTsolve");
401  sc = v_resize(sc,A->m);
402  n = sc->dim;
403  p = c->dim;
404  if ( n == p )
405  k = p-2;
406  else
407  k = p-1;
408  v_zero(sc);
409  sc->ve[0] = c->ve[0]/A->me[0][0];
410  if ( n == 1)
411  return sc;
412  if ( p > 1)
413  {
414  for ( i = 1; i < p; i++ )
415  {
416  s = 0.0;
417  for ( j = 0; j < i; j++ )
418  s += A->me[j][i]*sc->ve[j];
419  if ( A->me[i][i] == 0.0 )
420  error(E_SING,"QRTsolve");
421  sc->ve[i]=(c->ve[i]-s)/A->me[i][i];
422  }
423  }
424  for (i = k; i >= 0; i--)
425  {
426  s = diag->ve[i]*sc->ve[i];
427  for ( j = i+1; j < n; j++ )
428  s += A->me[j][i]*sc->ve[j];
429  r_ii = fabs(A->me[i][i]);
430  tmp_val = (r_ii*fabs(diag->ve[i]));
431  beta = ( tmp_val == 0.0 ) ? 0.0 : 1.0/tmp_val;
432  tmp_val = beta*s;
433  sc->ve[i] -= tmp_val*diag->ve[i];
434  for ( j = i+1; j < n; j++ )
435  sc->ve[j] -= tmp_val*A->me[j][i];
436  }
437 
438  return sc;
439 }
440 
441 /* QRcondest -- returns an estimate of the 2-norm condition number of the
442  matrix factorised by QRfactor() or QRCPfactor()
443  -- note that as Q does not affect the 2-norm condition number,
444  it is not necessary to pass the diag, beta (or pivot) vectors
445  -- generates a lower bound on the true condition number
446  -- if the matrix is exactly singular, HUGE is returned
447  -- note that QRcondest() is likely to be more reliable for
448  matrices factored using QRCPfactor() */
449 double QRcondest(QR)
450 MAT *QR;
451 {
452  static VEC *y=VNULL;
453  Real norm1, norm2, sum, tmp1, tmp2;
454  int i, j, limit;
455 
456  if ( QR == MNULL )
457  error(E_NULL,"QRcondest");
458 
459  limit = min(QR->m,QR->n);
460  for ( i = 0; i < limit; i++ )
461  if ( QR->me[i][i] == 0.0 )
462  return HUGE;
463 
464  y = v_resize(y,limit);
466  /* use the trick for getting a unit vector y with ||R.y||_inf small
467  from the LU condition estimator */
468  for ( i = 0; i < limit; i++ )
469  {
470  sum = 0.0;
471  for ( j = 0; j < i; j++ )
472  sum -= QR->me[j][i]*y->ve[j];
473  sum -= (sum < 0.0) ? 1.0 : -1.0;
474  y->ve[i] = sum / QR->me[i][i];
475  }
476  UTmlt(QR,y,y);
477 
478  /* now apply inverse power method to R^T.R */
479  for ( i = 0; i < 3; i++ )
480  {
481  tmp1 = v_norm2(y);
482  sv_mlt(1/tmp1,y,y);
483  UTsolve(QR,y,y,0.0);
484  tmp2 = v_norm2(y);
485  sv_mlt(1/v_norm2(y),y,y);
486  Usolve(QR,y,y,0.0);
487  }
488  /* now compute approximation for ||R^{-1}||_2 */
489  norm1 = sqrt(tmp1)*sqrt(tmp2);
490 
491  /* now use complementary approach to compute approximation to ||R||_2 */
492  for ( i = limit-1; i >= 0; i-- )
493  {
494  sum = 0.0;
495  for ( j = i+1; j < limit; j++ )
496  sum += QR->me[i][j]*y->ve[j];
497  y->ve[i] = (sum >= 0.0) ? 1.0 : -1.0;
498  y->ve[i] = (QR->me[i][i] >= 0.0) ? y->ve[i] : - y->ve[i];
499  }
500 
501  /* now apply power method to R^T.R */
502  for ( i = 0; i < 3; i++ )
503  {
504  tmp1 = v_norm2(y);
505  sv_mlt(1/tmp1,y,y);
506  Umlt(QR,y,y);
507  tmp2 = v_norm2(y);
508  sv_mlt(1/tmp2,y,y);
509  UTmlt(QR,y,y);
510  }
511  norm2 = sqrt(tmp1)*sqrt(tmp2);
512 
513  /* printf("QRcondest: norm1 = %g, norm2 = %g\n",norm1,norm2); */
514 
515  return norm1*norm2;
516 }
#define c
#define TRUE
Definition: err.c:57
#define FALSE
Definition: err.c:56
#define error(err_num, fn_name)
Definition: err.h:73
#define E_NULL
Definition: err.h:102
#define E_SING
Definition: err.h:98
#define E_SIZES
Definition: err.h:95
#define diag(s)
Definition: fmenu.cpp:192
MAT * hhtrcols(MAT *M, u_int i0, u_int j0, VEC *hh, double beta)
Definition: hsehldr.c:146
VEC * hhtrvec(VEC *hh, double beta, u_int i0, VEC *in, VEC *out)
Definition: hsehldr.c:72
VEC * hhvec(VEC *vec, u_int i0, Real *beta, VEC *out, Real *newval)
Definition: hsehldr.c:47
VEC * v_zero(VEC *x)
Definition: init.c:40
static Object ** v_resize(void *v)
Definition: ivocvect.cpp:1383
double __ip__(Real *dp1, Real *dp2, int len)
Definition: machine.c:40
#define Real
Definition: machine.h:189
uint32_t u_int
Definition: machine.h:38
VEC * UTsolve(MAT *A, VEC *b, VEC *x, double diag_val)
PERM * px_transp(PERM *px, u_int i, u_int j)
Definition: pxop.c:208
#define VNULL
Definition: matrix.h:631
#define v_copy(in, out)
Definition: matrix.h:404
#define m_copy(in, out)
Definition: matrix.h:403
#define set_col(mat, col, vec)
Definition: matrix.h:564
#define min(a, b)
Definition: matrix.h:157
#define MNULL
Definition: matrix.h:632
VEC * pxinv_vec(PERM *, VEC *, VEC *)
VEC * sv_mlt(double, VEC *, VEC *)
VEC * get_col(MAT *, u_int, VEC *)
#define v_norm2(x)
Definition: matrix.h:511
#define V_FREE(vec)
Definition: matrix.h:214
#define i
Definition: md1redef.h:12
#define TYPE_VEC
Definition: meminfo.h:52
#define MEM_STAT_REG(var, type)
Definition: meminfo.h:133
sqrt
Definition: extdef.h:3
fabs
Definition: extdef.h:3
#define A(i)
Definition: multisplit.cpp:63
int const size_t const size_t n
Definition: nrngsl.h:11
for(i=0;i< n;i++)
if(status)
size_t p
size_t j
static philox4x32_key_t k
Definition: nrnran123.cpp:11
MAT * m_get(int, int)
Definition: memory.c:36
VEC * Usolve()
MAT * makeQ(MAT *QR, diag, MAT *Qout)
Definition: qrfactor.c:226
MAT * makeR(MAT *QR, MAT *Rout)
Definition: qrfactor.c:276
static VEC * Umlt(MAT *U, VEC *x, VEC *out)
Definition: qrfactor.c:343
MAT * QRfactor(MAT *A, diag)
Definition: qrfactor.c:65
VEC * _Qsolve(MAT *QR, diag, VEC *diag *b, VEC *diag *x, VEC *diag *tmp)
Definition: qrfactor.c:185
static VEC * UTmlt(MAT *U, VEC *x, VEC *out)
Definition: qrfactor.c:363
double QRcondest(MAT *QR)
Definition: qrfactor.c:449
MAT * QRCPfactor(MAT *A, diag, PERM *px)
Definition: qrfactor.c:102
VEC * QRTsolve(MAT *A, VEC *diag, VEC *c, VEC *sc)
Definition: qrfactor.c:390
VEC * QRCPsolve(MAT *QR, diag, PERM *pivot, VEC *b, VEC *x)
Definition: qrfactor.c:321
VEC * QRsolve(MAT *QR, diag, VEC *diag *b, VEC *diag *x)
Definition: qrfactor.c:295
static char rcsid[]
Definition: qrfactor.c:42
#define NULL
Definition: sptree.h:16
Definition: matrix.h:73
Definition: matrix.h:87
u_int * pe
Definition: matrix.h:88
u_int size
Definition: matrix.h:88
Definition: matrix.h:67
u_int dim
Definition: matrix.h:68
Real * ve
Definition: matrix.h:69
#define square(x)
Definition: znorm.c:75