NEURON
itersym.c
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 
3 /**************************************************************************
4 **
5 ** Copyright (C) 1993 David E. Stewart & 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 /* itersym.c 17/09/93 */
29 
30 
31 /*
32  ITERATIVE METHODS - implementation of several iterative methods;
33  see also iter0.c
34  */
35 
36 #include <stdio.h>
37 #include "matrix.h"
38 #include "matrix2.h"
39 #include "sparse.h"
40 #include "iter.h"
41 #include <math.h>
42 
43 static char rcsid[] = "itersym.c,v 1.1 1997/12/04 17:55:29 hines Exp";
44 
45 
46 #ifdef ANSI_C
47 VEC *spCHsolve(SPMAT *,VEC *,VEC *);
48 VEC *trieig(VEC *,VEC *,MAT *);
49 #else
50 VEC *spCHsolve();
51 VEC *trieig();
52 #endif
53 
54 
55 
56 /* iter_spcg -- a simple interface to iter_cg() which uses sparse matrix
57  data structures
58  -- assumes that LLT contains the Cholesky factorisation of the
59  actual preconditioner;
60  use always as follows:
61  x = iter_spcg(A,LLT,b,eps,x,limit,steps);
62  or
63  x = iter_spcg(A,LLT,b,eps,VNULL,limit,steps);
64  In the second case the solution vector is created.
65  */
67 SPMAT *A, *LLT;
68 VEC *b, *x;
69 double eps;
70 int *steps, limit;
71 {
72  ITER *ip;
73 
74  ip = iter_get(0,0);
75  ip->Ax = (Fun_Ax) sp_mv_mlt;
76  ip->A_par = (void *)A;
77  ip->Bx = (Fun_Ax) spCHsolve;
78  ip->B_par = (void *)LLT;
79  ip->info = (Fun_info) NULL;
80  ip->b = b;
81  ip->eps = eps;
82  ip->limit = limit;
83  ip->x = x;
84  iter_cg(ip);
85  x = ip->x;
86  if (steps) *steps = ip->steps;
87  ip->shared_x = ip->shared_b = TRUE;
88  iter_free(ip); /* release only ITER structure */
89  return x;
90 }
91 
92 /*
93  Conjugate gradients method;
94  */
96 ITER *ip;
97 {
98  static VEC *r = VNULL, *p = VNULL, *q = VNULL, *z = VNULL;
99  Real alpha, beta, inner, old_inner, nres;
100  VEC *rr; /* rr == r or rr == z */
101 
102  if (ip == INULL)
103  error(E_NULL,"iter_cg");
104  if (!ip->Ax || !ip->b)
105  error(E_NULL,"iter_cg");
106  if ( ip->x == ip->b )
107  error(E_INSITU,"iter_cg");
108  if (!ip->stop_crit)
109  error(E_NULL,"iter_cg");
110 
111  if ( ip->eps <= 0.0 )
112  ip->eps = MACHEPS;
113 
114  r = v_resize(r,ip->b->dim);
115  p = v_resize(p,ip->b->dim);
116  q = v_resize(q,ip->b->dim);
117 
121 
122  if (ip->Bx != (Fun_Ax)NULL) {
123  z = v_resize(z,ip->b->dim);
125  rr = z;
126  }
127  else rr = r;
128 
129  if (ip->x != VNULL) {
130  if (ip->x->dim != ip->b->dim)
131  error(E_SIZES,"iter_cg");
132  ip->Ax(ip->A_par,ip->x,p); /* p = A*x */
133  v_sub(ip->b,p,r); /* r = b - A*x */
134  }
135  else { /* ip->x == 0 */
136  ip->x = v_get(ip->b->dim);
137  ip->shared_x = FALSE;
138  v_copy(ip->b,r);
139  }
140 
141  old_inner = 0.0;
142  for ( ip->steps = 0; ip->steps <= ip->limit; ip->steps++ )
143  {
144  if ( ip->Bx )
145  (ip->Bx)(ip->B_par,r,rr); /* rr = B*r */
146 
147  inner = in_prod(rr,r);
148  nres = sqrt(fabs(inner));
149  if (ip->info) ip->info(ip,nres,r,rr);
150  if (ip->steps == 0) ip->init_res = nres;
151  if ( ip->stop_crit(ip,nres,r,rr) ) break;
152 
153  if ( ip->steps ) /* if ( ip->steps > 0 ) ... */
154  {
155  beta = inner/old_inner;
156  p = v_mltadd(rr,p,beta,p);
157  }
158  else /* if ( ip->steps == 0 ) ... */
159  {
160  beta = 0.0;
161  p = v_copy(rr,p);
162  old_inner = 0.0;
163  }
164  (ip->Ax)(ip->A_par,p,q); /* q = A*p */
165  alpha = in_prod(p,q);
166  if (sqrt(fabs(alpha)) <= MACHEPS*ip->init_res)
167  error(E_BREAKDOWN,"iter_cg");
168  alpha = inner/alpha;
169  v_mltadd(ip->x,p,alpha,ip->x);
170  v_mltadd(r,q,-alpha,r);
171  old_inner = inner;
172  }
173 
174  return ip->x;
175 }
176 
177 
178 
179 /* iter_lanczos -- raw lanczos algorithm -- no re-orthogonalisation
180  -- creates T matrix of size == m,
181  but no larger than before beta_k == 0
182  -- uses passed routine to do matrix-vector multiplies */
183 void iter_lanczos(ip,a,b,beta2,Q)
184 ITER *ip;
185 VEC *a, *b;
186 Real *beta2;
187 MAT *Q;
188 {
189  int j;
190  static VEC *v = VNULL, *w = VNULL, *tmp = VNULL;
191  Real alpha, beta, c;
192 
193  if ( ! ip )
194  error(E_NULL,"iter_lanczos");
195  if ( ! ip->Ax || ! ip->x || ! a || ! b )
196  error(E_NULL,"iter_lanczos");
197  if ( ip->k <= 0 )
198  error(E_BOUNDS,"iter_lanczos");
199  if ( Q && ( Q->n < ip->x->dim || Q->m < ip->k ) )
200  error(E_SIZES,"iter_lanczos");
201 
202  a = v_resize(a,(u_int)ip->k);
203  b = v_resize(b,(u_int)(ip->k-1));
204  v = v_resize(v,ip->x->dim);
205  w = v_resize(w,ip->x->dim);
206  tmp = v_resize(tmp,ip->x->dim);
209  MEM_STAT_REG(tmp,TYPE_VEC);
210 
211  beta = 1.0;
212  v_zero(a);
213  v_zero(b);
214  if (Q) m_zero(Q);
215 
216  /* normalise x as w */
217  c = v_norm2(ip->x);
218  if (c <= MACHEPS) { /* ip->x == 0 */
219  *beta2 = 0.0;
220  return;
221  }
222  else
223  sv_mlt(1.0/c,ip->x,w);
224 
225  (ip->Ax)(ip->A_par,w,v);
226 
227  for ( j = 0; j < ip->k; j++ )
228  {
229  /* store w in Q if Q not NULL */
230  if ( Q ) set_row(Q,j,w);
231 
232  alpha = in_prod(w,v);
233  a->ve[j] = alpha;
234  v_mltadd(v,w,-alpha,v);
235  beta = v_norm2(v);
236  if ( beta == 0.0 )
237  {
238  *beta2 = 0.0;
239  return;
240  }
241 
242  if ( j < ip->k-1 )
243  b->ve[j] = beta;
244  v_copy(w,tmp);
245  sv_mlt(1/beta,v,w);
246  sv_mlt(-beta,tmp,v);
247  (ip->Ax)(ip->A_par,w,tmp);
248  v_add(v,tmp,v);
249  }
250  *beta2 = beta;
251 
252 }
253 
254 /* iter_splanczos -- version that uses sparse matrix data structure */
255 void iter_splanczos(A,m,x0,a,b,beta2,Q)
256 SPMAT *A;
257 int m;
258 VEC *x0, *a, *b;
259 Real *beta2;
260 MAT *Q;
261 {
262  ITER *ip;
263 
264  ip = iter_get(0,0);
265  ip->shared_x = ip->shared_b = TRUE;
266  ip->Ax = (Fun_Ax) sp_mv_mlt;
267  ip->A_par = (void *) A;
268  ip->x = x0;
269  ip->k = m;
270  iter_lanczos(ip,a,b,beta2,Q);
271  iter_free(ip); /* release only ITER structure */
272 }
273 
274 
275 #ifndef MAC
276 extern double frexp(), ldexp();
277 #endif
278 
279 /* product -- returns the product of a long list of numbers
280  -- answer stored in mant (mantissa) and expt (exponent) */
281 static double product(a,offset,expt)
282 VEC *a;
283 double offset;
284 int *expt;
285 {
286  Real mant, tmp_fctr;
287  int i, tmp_expt;
288 
289  if ( ! a )
290  error(E_NULL,"product");
291 
292  mant = 1.0;
293  *expt = 0;
294  if ( offset == 0.0 )
295  for ( i = 0; i < a->dim; i++ )
296  {
297  mant *= frexp(a->ve[i],&tmp_expt);
298  *expt += tmp_expt;
299  if ( ! (i % 10) )
300  {
301  mant = frexp(mant,&tmp_expt);
302  *expt += tmp_expt;
303  }
304  }
305  else
306  for ( i = 0; i < a->dim; i++ )
307  {
308  tmp_fctr = a->ve[i] - offset;
309  tmp_fctr += (tmp_fctr > 0.0 ) ? -MACHEPS*offset :
310  MACHEPS*offset;
311  mant *= frexp(tmp_fctr,&tmp_expt);
312  *expt += tmp_expt;
313  if ( ! (i % 10) )
314  {
315  mant = frexp(mant,&tmp_expt);
316  *expt += tmp_expt;
317  }
318  }
319 
320  mant = frexp(mant,&tmp_expt);
321  *expt += tmp_expt;
322 
323  return mant;
324 }
325 
326 /* product2 -- returns the product of a long list of numbers
327  -- answer stored in mant (mantissa) and expt (exponent) */
328 static double product2(a,k,expt)
329 VEC *a;
330 int k; /* entry of a to leave out */
331 int *expt;
332 {
333  Real mant, mu, tmp_fctr;
334  int i, tmp_expt;
335 
336  if ( ! a )
337  error(E_NULL,"product2");
338  if ( k < 0 || k >= a->dim )
339  error(E_BOUNDS,"product2");
340 
341  mant = 1.0;
342  *expt = 0;
343  mu = a->ve[k];
344  for ( i = 0; i < a->dim; i++ )
345  {
346  if ( i == k )
347  continue;
348  tmp_fctr = a->ve[i] - mu;
349  tmp_fctr += ( tmp_fctr > 0.0 ) ? -MACHEPS*mu : MACHEPS*mu;
350  mant *= frexp(tmp_fctr,&tmp_expt);
351  *expt += tmp_expt;
352  if ( ! (i % 10) )
353  {
354  mant = frexp(mant,&tmp_expt);
355  *expt += tmp_expt;
356  }
357  }
358  mant = frexp(mant,&tmp_expt);
359  *expt += tmp_expt;
360 
361  return mant;
362 }
363 
364 /* dbl_cmp -- comparison function to pass to qsort() */
365 static int dbl_cmp(x,y)
366 Real *x, *y;
367 {
368  Real tmp;
369 
370  tmp = *x - *y;
371  return (tmp > 0 ? 1 : tmp < 0 ? -1: 0);
372 }
373 
374 /* iter_lanczos2 -- lanczos + error estimate for every e-val
375  -- uses Cullum & Willoughby approach, Sparse Matrix Proc. 1978
376  -- returns multiple e-vals where multiple e-vals may not exist
377  -- returns evals vector */
378 VEC *iter_lanczos2(ip,evals,err_est)
379 ITER *ip; /* ITER structure */
380 VEC *evals; /* eigenvalue vector */
381 VEC *err_est; /* error estimates of eigenvalues */
382 {
383  VEC *a;
384  static VEC *b=VNULL, *a2=VNULL, *b2=VNULL;
385  Real beta, pb_mant, det_mant, det_mant1, det_mant2;
386  int i, pb_expt, det_expt, det_expt1, det_expt2;
387 
388  if ( ! ip )
389  error(E_NULL,"iter_lanczos2");
390  if ( ! ip->Ax || ! ip->x )
391  error(E_NULL,"iter_lanczos2");
392  if ( ip->k <= 0 )
393  error(E_RANGE,"iter_lanczos2");
394 
395  a = evals;
396  a = v_resize(a,(u_int)ip->k);
397  b = v_resize(b,(u_int)(ip->k-1));
399 
400  iter_lanczos(ip,a,b,&beta,MNULL);
401 
402  /* printf("# beta =%g\n",beta); */
403  pb_mant = 0.0;
404  if ( err_est )
405  {
406  pb_mant = product(b,(double)0.0,&pb_expt);
407  /* printf("# pb_mant = %g, pb_expt = %d\n",pb_mant, pb_expt); */
408  }
409 
410  /* printf("# diags =\n"); v_output(a); */
411  /* printf("# off diags =\n"); v_output(b); */
412  a2 = v_resize(a2,a->dim - 1);
413  b2 = v_resize(b2,b->dim - 1);
416  for ( i = 0; i < a2->dim - 1; i++ )
417  {
418  a2->ve[i] = a->ve[i+1];
419  b2->ve[i] = b->ve[i+1];
420  }
421  a2->ve[a2->dim-1] = a->ve[a2->dim];
422 
423  trieig(a,b,MNULL);
424 
425  /* sort evals as a courtesy */
426  qsort((void *)(a->ve),(int)(a->dim),sizeof(Real),(int (*)())dbl_cmp);
427 
428  /* error estimates */
429  if ( err_est )
430  {
431  err_est = v_resize(err_est,(u_int)ip->k);
432 
433  trieig(a2,b2,MNULL);
434  /* printf("# a =\n"); v_output(a); */
435  /* printf("# a2 =\n"); v_output(a2); */
436 
437  for ( i = 0; i < a->dim; i++ )
438  {
439  det_mant1 = product2(a,i,&det_expt1);
440  det_mant2 = product(a2,(double)a->ve[i],&det_expt2);
441  /* printf("# det_mant1=%g, det_expt1=%d\n",
442  det_mant1,det_expt1); */
443  /* printf("# det_mant2=%g, det_expt2=%d\n",
444  det_mant2,det_expt2); */
445  if ( det_mant1 == 0.0 )
446  { /* multiple e-val of T */
447  err_est->ve[i] = 0.0;
448  continue;
449  }
450  else if ( det_mant2 == 0.0 )
451  {
452  err_est->ve[i] = HUGE;
453  continue;
454  }
455  if ( (det_expt1 + det_expt2) % 2 )
456  /* if odd... */
457  det_mant = sqrt(2.0*fabs(det_mant1*det_mant2));
458  else /* if even... */
459  det_mant = sqrt(fabs(det_mant1*det_mant2));
460  det_expt = (det_expt1+det_expt2)/2;
461  err_est->ve[i] = fabs(beta*
462  ldexp(pb_mant/det_mant,pb_expt-det_expt));
463  }
464  }
465 
466  return a;
467 }
468 
469 /* iter_splanczos2 -- version of iter_lanczos2() that uses sparse matrix data
470  structure */
471 
472 VEC *iter_splanczos2(A,m,x0,evals,err_est)
473 SPMAT *A;
474 int m;
475 VEC *x0; /* initial vector */
476 VEC *evals; /* eigenvalue vector */
477 VEC *err_est; /* error estimates of eigenvalues */
478 {
479  ITER *ip;
480  VEC *a;
481 
482  ip = iter_get(0,0);
483  ip->Ax = (Fun_Ax) sp_mv_mlt;
484  ip->A_par = (void *) A;
485  ip->x = x0;
486  ip->k = m;
487  a = iter_lanczos2(ip,evals,err_est);
488  ip->shared_x = ip->shared_b = TRUE;
489  iter_free(ip); /* release only ITER structure */
490  return a;
491 }
492 
493 
494 
495 
496 /*
497  Conjugate gradient method
498  Another variant - mainly for testing
499  */
500 
502 ITER *ip;
503 {
504  static VEC *r = VNULL, *p = VNULL, *q = VNULL, *z = VNULL;
505  Real alpha;
506  double inner,nres;
507  VEC *rr; /* rr == r or rr == z */
508 
509  if (ip == INULL)
510  error(E_NULL,"iter_cg");
511  if (!ip->Ax || !ip->b)
512  error(E_NULL,"iter_cg");
513  if ( ip->x == ip->b )
514  error(E_INSITU,"iter_cg");
515  if (!ip->stop_crit)
516  error(E_NULL,"iter_cg");
517 
518  if ( ip->eps <= 0.0 )
519  ip->eps = MACHEPS;
520 
521  r = v_resize(r,ip->b->dim);
522  p = v_resize(p,ip->b->dim);
523  q = v_resize(q,ip->b->dim);
524 
528 
529  if (ip->Bx != (Fun_Ax)NULL) {
530  z = v_resize(z,ip->b->dim);
532  rr = z;
533  }
534  else rr = r;
535 
536  if (ip->x != VNULL) {
537  if (ip->x->dim != ip->b->dim)
538  error(E_SIZES,"iter_cg");
539  ip->Ax(ip->A_par,ip->x,p); /* p = A*x */
540  v_sub(ip->b,p,r); /* r = b - A*x */
541  }
542  else { /* ip->x == 0 */
543  ip->x = v_get(ip->b->dim);
544  ip->shared_x = FALSE;
545  v_copy(ip->b,r);
546  }
547 
548  if (ip->Bx) (ip->Bx)(ip->B_par,r,p);
549  else v_copy(r,p);
550 
551  inner = in_prod(p,r);
552  nres = sqrt(fabs(inner));
553  if (ip->info) ip->info(ip,nres,r,p);
554  if ( nres == 0.0) return ip->x;
555 
556  for ( ip->steps = 0; ip->steps <= ip->limit; ip->steps++ )
557  {
558  ip->Ax(ip->A_par,p,q);
559  inner = in_prod(q,p);
560  if (sqrt(fabs(inner)) <= MACHEPS*ip->init_res)
561  error(E_BREAKDOWN,"iter_cg1");
562 
563  alpha = in_prod(p,r)/inner;
564  v_mltadd(ip->x,p,alpha,ip->x);
565  v_mltadd(r,q,-alpha,r);
566 
567  rr = r;
568  if (ip->Bx) {
569  ip->Bx(ip->B_par,r,z);
570  rr = z;
571  }
572 
573  nres = in_prod(r,rr);
574  if (nres < 0.0) {
575  warning(WARN_RES_LESS_0,"iter_cg");
576  break;
577  }
578  nres = sqrt(fabs(nres));
579  if (ip->info) ip->info(ip,nres,r,z);
580  if (ip->steps == 0) ip->init_res = nres;
581  if ( ip->stop_crit(ip,nres,r,z) ) break;
582 
583  alpha = -in_prod(rr,q)/inner;
584  v_mltadd(rr,p,alpha,p);
585 
586  }
587 
588  return ip->x;
589 }
590 
591 
Fun_Ax Ax
Definition: iter.h:68
void iter_splanczos(SPMAT *A, int m, VEC *x0, VEC *a, VEC *b, Real *beta2, MAT *Q)
Definition: itersym.c:255
void(* info)()
Definition: iter.h:92
#define WARN_RES_LESS_0
Definition: err.h:122
VEC * iter_splanczos2(SPMAT *A, int m, VEC *x0, VEC *evals, VEC *err_est)
Definition: itersym.c:472
unsigned k
Definition: iter.h:59
#define TYPE_VEC
Definition: meminfo.h:52
#define in_prod(a, b)
Definition: matrix.h:485
VEC * iter_spcg(SPMAT *A, SPMAT *LLT, VEC *b, double eps, VEC *x, int limit, int *steps)
Definition: itersym.c:66
#define Real
Definition: machine.h:189
static Object ** v_resize(void *v)
Definition: ivocvect.cpp:1297
double ldexp()
size_t p
VEC * spCHsolve(SPMAT *, VEC *, VEC *)
Definition: spchfctr.c:311
#define TRUE
Definition: err.c:57
#define v_norm2(x)
Definition: matrix.h:511
VEC * trieig(VEC *, VEC *, MAT *)
Definition: symmeig.c:49
#define v
Definition: md1redef.h:4
int iter_free(ITER *ip)
Definition: iter0.c:118
static philox4x32_key_t k
Definition: nrnran123.cpp:11
u_int n
Definition: matrix.h:74
int shared_b
Definition: iter.h:58
int steps
Definition: iter.h:61
ITER * iter_get(int lenb, int lenx)
Definition: iter0.c:76
Real * ve
Definition: matrix.h:69
Definition: matrix.h:67
void iter_lanczos(ITER *ip, VEC *a, VEC *b, Real *beta2, MAT *Q)
Definition: itersym.c:183
static Object ** v_sub(void *v1)
Definition: ivocvect.cpp:2263
VEC * iter_lanczos2(ITER *ip, VEC *evals, VEC *err_est)
Definition: itersym.c:378
#define v_copy(in, out)
Definition: matrix.h:404
static double product(VEC *a, double offset, int *expt)
Definition: itersym.c:281
#define E_SIZES
Definition: err.h:95
uint32_t u_int
Definition: machine.h:38
Definition: sparse.h:54
#define E_RANGE
Definition: err.h:104
VEC * iter_cg(ITER *ip)
Definition: itersym.c:95
static char rcsid[]
Definition: itersym.c:43
#define MNULL
Definition: matrix.h:632
u_int dim
Definition: matrix.h:68
#define INULL
Definition: iter.h:113
int
Definition: nrnmusic.cpp:71
Fun_Ax Bx
Definition: iter.h:75
VEC *(* Fun_Ax)(void *, VEC *, VEC *)
Definition: iter.h:49
VEC * iter_cg1(ITER *ip)
Definition: itersym.c:501
double frexp()
#define E_BOUNDS
Definition: err.h:96
VEC * v_mltadd(VEC *, VEC *, double, VEC *)
#define E_NULL
Definition: err.h:102
size_t j
#define alpha
Definition: bkpfacto.c:43
int limit
Definition: iter.h:60
VEC * sv_mlt(double, VEC *, VEC *)
VEC * sp_mv_mlt(SPMAT *A, VEC *x, VEC *out)
Definition: sparse.c:123
static double v_get(void *v)
Definition: ivocvect.cpp:1309
static double product2(VEC *a, int k, int *expt)
Definition: itersym.c:328
#define MACHEPS
Definition: machine.h:213
VEC * v_zero(VEC *x)
Definition: init.c:40
int shared_x
Definition: iter.h:57
void warning(const char *s, const char *t)
Definition: hoc.cpp:1542
sqrt
Definition: extdef.h:3
Real eps
Definition: iter.h:62
static Object ** v_add(void *v1)
Definition: ivocvect.cpp:2245
#define E_BREAKDOWN
Definition: err.h:116
#define MEM_STAT_REG(var, type)
Definition: meminfo.h:133
#define FALSE
Definition: err.c:56
#define i
Definition: md1redef.h:12
#define A(i)
Definition: multisplit.cpp:61
#define c
static Object ** m_zero(void *v)
Definition: matrix.cpp:511
#define error(err_num, fn_name)
Definition: err.h:73
fabs
Definition: extdef.h:3
#define set_row(mat, row, vec)
Definition: matrix.h:562
u_int m
Definition: matrix.h:74
void(* Fun_info)(ITER *, double, VEC *, VEC *)
Definition: iter.h:117
VEC * x
Definition: iter.h:64
#define VNULL
Definition: matrix.h:631
static int dbl_cmp(Real *x, Real *y)
Definition: itersym.c:365
size_t q
void * B_par
Definition: iter.h:76
return NULL
Definition: cabcode.cpp:461
Definition: iter.h:56
Definition: matrix.h:73
void * A_par
Definition: iter.h:69
VEC * b
Definition: iter.h:66
#define E_INSITU
Definition: err.h:106