NEURON
spbkp.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  Sparse matrix Bunch--Kaufman--Parlett factorisation and solve
30  Radical revision started Thu 05th Nov 1992, 09:36:12 AM
31  to use Karen George's suggestion of leaving the the row elements unordered
32  Radical revision completed Mon 07th Dec 1992, 10:59:57 AM
33 */
34 
35 static char rcsid[] = "spbkp.c,v 1.1 1997/12/04 17:55:50 hines Exp";
36 
37 #include <stdio.h>
38 #include "sparse2.h"
39 #include <math.h>
40 
41 
42 #ifdef MALLOCDECL
43 #include <malloc.h>
44 #endif
45 
46 #define alpha 0.6403882032022076 /* = (1+sqrt(17))/8 */
47 
48 
49 #define btos(x) ((x) ? "TRUE" : "FALSE")
50 
51 /* assume no use of sqr() uses side-effects */
52 #define sqr(x) ((x)*(x))
53 
54 /* unord_get_idx -- returns index (encoded if entry not allocated)
55  of the element of row r with column j
56  -- uses linear search */
58 SPROW *r;
59 int j;
60 {
61  int idx;
62  row_elt *e;
63 
64  if ( ! r || ! r->elt )
65  error(E_NULL,"unord_get_idx");
66  for ( idx = 0, e = r->elt; idx < r->len; idx++, e++ )
67  if ( e->col == j )
68  break;
69  if ( idx >= r->len )
70  return -(r->len+2);
71  else
72  return idx;
73 }
74 
75 /* unord_get_val -- returns value of the (i,j) entry of A
76  -- same assumptions as unord_get_idx() */
77 double unord_get_val(A,i,j)
78 SPMAT *A;
79 int i, j;
80 {
81  SPROW *r;
82  int idx;
83 
84  if ( ! A )
85  error(E_NULL,"unord_get_val");
86  if ( i < 0 || i >= A->m || j < 0 || j >= A->n )
87  error(E_BOUNDS,"unord_get_val");
88 
89  r = &(A->row[i]);
90  idx = unord_get_idx(r,j);
91  if ( idx < 0 )
92  return 0.0;
93  else
94  return r->elt[idx].val;
95 }
96 
97 
98 /* bkp_swap_elt -- swaps the (i,j) with the (k,l) entry of sparse matrix
99  -- either or both of the entries may be unallocated */
100 static SPMAT *bkp_swap_elt(A,i1,j1,idx1,i2,j2,idx2)
101 SPMAT *A;
102 int i1, j1, idx1, i2, j2, idx2;
103 {
104  int tmp_row, tmp_idx;
105  SPROW *r1, *r2;
106  row_elt *e1, *e2;
107  Real tmp;
108 
109  if ( ! A )
110  error(E_NULL,"bkp_swap_elt");
111 
112  if ( i1 < 0 || j1 < 0 || i2 < 0 || j2 < 0 ||
113  i1 >= A->m || j1 >= A->n || i2 >= A->m || j2 >= A->n )
114  {
115  error(E_BOUNDS,"bkp_swap_elt");
116  }
117 
118  if ( i1 == i2 && j1 == j2 )
119  return A;
120  if ( idx1 < 0 && idx2 < 0 ) /* neither allocated */
121  return A;
122 
123  r1 = &(A->row[i1]); r2 = &(A->row[i2]);
124  /* if ( idx1 >= r1->len || idx2 >= r2->len )
125  error(E_BOUNDS,"bkp_swap_elt"); */
126  if ( idx1 < 0 ) /* assume not allocated */
127  {
128  idx1 = r1->len;
129  if ( idx1 >= r1->maxlen )
130  { tracecatch(sprow_xpd(r1,2*r1->maxlen+1,TYPE_SPMAT),
131  "bkp_swap_elt"); }
132  r1->len = idx1+1;
133  r1->elt[idx1].col = j1;
134  r1->elt[idx1].val = 0.0;
135  /* now patch up column access path */
136  tmp_row = -1; tmp_idx = j1;
137  chase_col(A,j1,&tmp_row,&tmp_idx,i1-1);
138 
139  if ( tmp_row < 0 )
140  {
141  r1->elt[idx1].nxt_row = A->start_row[j1];
142  r1->elt[idx1].nxt_idx = A->start_idx[j1];
143  A->start_row[j1] = i1;
144  A->start_idx[j1] = idx1;
145  }
146  else
147  {
148  row_elt *tmp_e;
149 
150  tmp_e = &(A->row[tmp_row].elt[tmp_idx]);
151  r1->elt[idx1].nxt_row = tmp_e->nxt_row;
152  r1->elt[idx1].nxt_idx = tmp_e->nxt_idx;
153  tmp_e->nxt_row = i1;
154  tmp_e->nxt_idx = idx1;
155  }
156  }
157  else if ( r1->elt[idx1].col != j1 )
158  error(E_INTERN,"bkp_swap_elt");
159  if ( idx2 < 0 )
160  {
161  idx2 = r2->len;
162  if ( idx2 >= r2->maxlen )
163  { tracecatch(sprow_xpd(r2,2*r2->maxlen+1,TYPE_SPMAT),
164  "bkp_swap_elt"); }
165 
166  r2->len = idx2+1;
167  r2->elt[idx2].col = j2;
168  r2->elt[idx2].val = 0.0;
169  /* now patch up column access path */
170  tmp_row = -1; tmp_idx = j2;
171  chase_col(A,j2,&tmp_row,&tmp_idx,i2-1);
172  if ( tmp_row < 0 )
173  {
174  r2->elt[idx2].nxt_row = A->start_row[j2];
175  r2->elt[idx2].nxt_idx = A->start_idx[j2];
176  A->start_row[j2] = i2;
177  A->start_idx[j2] = idx2;
178  }
179  else
180  {
181  row_elt *tmp_e;
182 
183  tmp_e = &(A->row[tmp_row].elt[tmp_idx]);
184  r2->elt[idx2].nxt_row = tmp_e->nxt_row;
185  r2->elt[idx2].nxt_idx = tmp_e->nxt_idx;
186  tmp_e->nxt_row = i2;
187  tmp_e->nxt_idx = idx2;
188  }
189  }
190  else if ( r2->elt[idx2].col != j2 )
191  error(E_INTERN,"bkp_swap_elt");
192 
193  e1 = &(r1->elt[idx1]); e2 = &(r2->elt[idx2]);
194 
195  tmp = e1->val;
196  e1->val = e2->val;
197  e2->val = tmp;
198 
199  return A;
200 }
201 
202 /* bkp_bump_col -- bumps row and idx to next entry in column j */
204 SPMAT *A;
205 int j, *row, *idx;
206 {
207  SPROW *r;
208  row_elt *e;
209 
210  if ( *row < 0 )
211  {
212  *row = A->start_row[j];
213  *idx = A->start_idx[j];
214  }
215  else
216  {
217  r = &(A->row[*row]);
218  e = &(r->elt[*idx]);
219  if ( e->col != j )
220  error(E_INTERN,"bkp_bump_col");
221  *row = e->nxt_row;
222  *idx = e->nxt_idx;
223  }
224  if ( *row < 0 )
225  return (row_elt *)NULL;
226  else
227  return &(A->row[*row].elt[*idx]);
228 }
229 
230 /* bkp_interchange -- swap rows/cols i and j (symmetric pivot)
231  -- uses just the upper triangular part */
233 SPMAT *A;
234 int i1, i2;
235 {
236  int tmp_row, tmp_idx;
237  int row1, row2, idx1, idx2, tmp_row1, tmp_idx1, tmp_row2, tmp_idx2;
238  SPROW *r1, *r2;
239  row_elt *e1, *e2;
240  IVEC *done_list = IVNULL;
241 
242  if ( ! A )
243  error(E_NULL,"bkp_interchange");
244  if ( i1 < 0 || i1 >= A->n || i2 < 0 || i2 >= A->n )
245  error(E_BOUNDS,"bkp_interchange");
246  if ( A->m != A->n )
247  error(E_SQUARE,"bkp_interchange");
248 
249  if ( i1 == i2 )
250  return A;
251  if ( i1 > i2 )
252  { tmp_idx = i1; i1 = i2; i2 = tmp_idx; }
253 
255  for ( tmp_idx = 0; tmp_idx < A->n; tmp_idx++ )
256  done_list->ive[tmp_idx] = FALSE;
257  row1 = -1; idx1 = i1;
258  row2 = -1; idx2 = i2;
259  e1 = bkp_bump_col(A,i1,&row1,&idx1);
260  e2 = bkp_bump_col(A,i2,&row2,&idx2);
261 
262  while ( (row1 >= 0 && row1 < i1) || (row2 >= 0 && row2 < i1) )
263  /* Note: "row2 < i1" not "row2 < i2" as we must stop before the
264  "knee bend" */
265  {
266  if ( row1 >= 0 && row1 < i1 && ( row1 < row2 || row2 < 0 ) )
267  {
268  tmp_row1 = row1; tmp_idx1 = idx1;
269  e1 = bkp_bump_col(A,i1,&tmp_row1,&tmp_idx1);
270  if ( ! done_list->ive[row1] )
271  {
272  if ( row1 == row2 )
273  bkp_swap_elt(A,row1,i1,idx1,row1,i2,idx2);
274  else
275  bkp_swap_elt(A,row1,i1,idx1,row1,i2,-1);
276  done_list->ive[row1] = TRUE;
277  }
278  row1 = tmp_row1; idx1 = tmp_idx1;
279  }
280  else if ( row2 >= 0 && row2 < i1 && ( row2 < row1 || row1 < 0 ) )
281  {
282  tmp_row2 = row2; tmp_idx2 = idx2;
283  e2 = bkp_bump_col(A,i2,&tmp_row2,&tmp_idx2);
284  if ( ! done_list->ive[row2] )
285  {
286  if ( row1 == row2 )
287  bkp_swap_elt(A,row2,i1,idx1,row2,i2,idx2);
288  else
289  bkp_swap_elt(A,row2,i1,-1,row2,i2,idx2);
290  done_list->ive[row2] = TRUE;
291  }
292  row2 = tmp_row2; idx2 = tmp_idx2;
293  }
294  else if ( row1 == row2 )
295  {
296  tmp_row1 = row1; tmp_idx1 = idx1;
297  e1 = bkp_bump_col(A,i1,&tmp_row1,&tmp_idx1);
298  tmp_row2 = row2; tmp_idx2 = idx2;
299  e2 = bkp_bump_col(A,i2,&tmp_row2,&tmp_idx2);
300  if ( ! done_list->ive[row1] )
301  {
302  bkp_swap_elt(A,row1,i1,idx1,row2,i2,idx2);
303  done_list->ive[row1] = TRUE;
304  }
305  row1 = tmp_row1; idx1 = tmp_idx1;
306  row2 = tmp_row2; idx2 = tmp_idx2;
307  }
308  }
309 
310  /* ensure we are **past** the first knee */
311  while ( row2 >= 0 && row2 <= i1 )
312  e2 = bkp_bump_col(A,i2,&row2,&idx2);
313 
314  /* at/after 1st "knee bend" */
315  r1 = &(A->row[i1]);
316  idx1 = 0;
317  e1 = &(r1->elt[idx1]);
318  while ( row2 >= 0 && row2 < i2 )
319  {
320  /* used for update of e2 at end of loop */
321  tmp_row = row2; tmp_idx = idx2;
322  if ( ! done_list->ive[row2] )
323  {
324  r2 = &(A->row[row2]);
325  bkp_bump_col(A,i2,&tmp_row,&tmp_idx);
326  done_list->ive[row2] = TRUE;
327  tmp_idx1 = unord_get_idx(r1,row2);
328  tracecatch(bkp_swap_elt(A,row2,i2,idx2,i1,row2,tmp_idx1),
329  "bkp_interchange");
330  }
331 
332  /* update e1 and e2 */
333  row2 = tmp_row; idx2 = tmp_idx;
334  e2 = ( row2 >= 0 ) ? &(A->row[row2].elt[idx2]) : (row_elt *)NULL;
335  }
336 
337  idx1 = 0;
338  e1 = r1->elt;
339  while ( idx1 < r1->len )
340  {
341  if ( e1->col >= i2 || e1->col <= i1 )
342  {
343  idx1++;
344  e1++;
345  continue;
346  }
347  if ( ! done_list->ive[e1->col] )
348  {
349  tmp_idx2 = unord_get_idx(&(A->row[e1->col]),i2);
350  tracecatch(bkp_swap_elt(A,i1,e1->col,idx1,e1->col,i2,tmp_idx2),
351  "bkp_interchange");
352  done_list->ive[e1->col] = TRUE;
353  }
354  idx1++;
355  e1++;
356  }
357 
358  /* at/after 2nd "knee bend" */
359  idx1 = 0;
360  e1 = &(r1->elt[idx1]);
361  r2 = &(A->row[i2]);
362  idx2 = 0;
363  e2 = &(r2->elt[idx2]);
364  while ( idx1 < r1->len )
365  {
366  if ( e1->col <= i2 )
367  {
368  idx1++; e1++;
369  continue;
370  }
371  if ( ! done_list->ive[e1->col] )
372  {
373  tmp_idx2 = unord_get_idx(r2,e1->col);
374  tracecatch(bkp_swap_elt(A,i1,e1->col,idx1,i2,e1->col,tmp_idx2),
375  "bkp_interchange");
376  done_list->ive[e1->col] = TRUE;
377  }
378  idx1++; e1++;
379  }
380 
381  idx2 = 0; e2 = r2->elt;
382  while ( idx2 < r2->len )
383  {
384  if ( e2->col <= i2 )
385  {
386  idx2++; e2++;
387  continue;
388  }
389  if ( ! done_list->ive[e2->col] )
390  {
391  tmp_idx1 = unord_get_idx(r1,e2->col);
392  tracecatch(bkp_swap_elt(A,i2,e2->col,idx2,i1,e2->col,tmp_idx1),
393  "bkp_interchange");
394  done_list->ive[e2->col] = TRUE;
395  }
396  idx2++; e2++;
397  }
398 
399  /* now interchange the digonal entries! */
400  idx1 = unord_get_idx(&(A->row[i1]),i1);
401  idx2 = unord_get_idx(&(A->row[i2]),i2);
402  if ( idx1 >= 0 || idx2 >= 0 )
403  {
404  tracecatch(bkp_swap_elt(A,i1,i1,idx1,i2,i2,idx2),
405  "bkp_interchange");
406  }
407 
408  return A;
409 }
410 
411 
412 /* iv_min -- returns minimum of an integer vector
413  -- sets index to the position in iv if index != NULL */
414 int iv_min(iv,index)
415 IVEC *iv;
416 int *index;
417 {
418  int i, i_min, min_val, tmp;
419 
420  if ( ! iv )
421  error(E_NULL,"iv_min");
422  if ( iv->dim <= 0 )
423  error(E_SIZES,"iv_min");
424  i_min = 0;
425  min_val = iv->ive[0];
426  for ( i = 1; i < iv->dim; i++ )
427  {
428  tmp = iv->ive[i];
429  if ( tmp < min_val )
430  {
431  min_val = tmp;
432  i_min = i;
433  }
434  }
435 
436  if ( index != (int *)NULL )
437  *index = i_min;
438 
439  return min_val;
440 }
441 
442 /* max_row_col -- returns max { |A[j][k]| : k >= i, k != j, k != l } given j
443  using symmetry and only the upper triangular part of A */
444 static double max_row_col(A,i,j,l)
445 SPMAT *A;
446 int i, j, l;
447 {
448  int row_num, idx;
449  SPROW *r;
450  row_elt *e;
451  Real max_val, tmp;
452 
453  if ( ! A )
454  error(E_NULL,"max_row_col");
455  if ( i < 0 || i > A->n || j < 0 || j >= A->n )
456  error(E_BOUNDS,"max_row_col");
457 
458  max_val = 0.0;
459 
460  idx = unord_get_idx(&(A->row[i]),j);
461  if ( idx < 0 )
462  {
463  row_num = -1; idx = j;
464  e = chase_past(A,j,&row_num,&idx,i);
465  }
466  else
467  {
468  row_num = i;
469  e = &(A->row[i].elt[idx]);
470  }
471  while ( row_num >= 0 && row_num < j )
472  {
473  if ( row_num != l )
474  {
475  tmp = fabs(e->val);
476  if ( tmp > max_val )
477  max_val = tmp;
478  }
479  e = bump_col(A,j,&row_num,&idx);
480  }
481  r = &(A->row[j]);
482  for ( idx = 0, e = r->elt; idx < r->len; idx++, e++ )
483  {
484  if ( e->col > j && e->col != l )
485  {
486  tmp = fabs(e->val);
487  if ( tmp > max_val )
488  max_val = tmp;
489  }
490  }
491 
492  return max_val;
493 }
494 
495 /* nonzeros -- counts non-zeros in A */
496 static int nonzeros(A)
497 SPMAT *A;
498 {
499  int cnt, i;
500 
501  if ( ! A )
502  return 0;
503  cnt = 0;
504  for ( i = 0; i < A->m; i++ )
505  cnt += A->row[i].len;
506 
507  return cnt;
508 }
509 
510 /* chk_col_access -- for spBKPfactor()
511  -- checks that column access path is OK */
513 SPMAT *A;
514 {
515  int cnt_nz, j, row, idx;
516  SPROW *r;
517  row_elt *e;
518 
519  if ( ! A )
520  error(E_NULL,"chk_col_access");
521 
522  /* count nonzeros as we go down columns */
523  cnt_nz = 0;
524  for ( j = 0; j < A->n; j++ )
525  {
526  row = A->start_row[j];
527  idx = A->start_idx[j];
528  while ( row >= 0 )
529  {
530  if ( row >= A->m || idx < 0 )
531  return FALSE;
532  r = &(A->row[row]);
533  if ( idx >= r->len )
534  return FALSE;
535  e = &(r->elt[idx]);
536  if ( e->nxt_row >= 0 && e->nxt_row <= row )
537  return FALSE;
538  row = e->nxt_row;
539  idx = e->nxt_idx;
540  cnt_nz++;
541  }
542  }
543 
544  if ( cnt_nz != nonzeros(A) )
545  return FALSE;
546  else
547  return TRUE;
548 }
549 
550 /* col_cmp -- compare two columns -- for sorting rows using qsort() */
551 static int col_cmp(e1,e2)
552 row_elt *e1, *e2;
553 {
554  return e1->col - e2->col;
555 }
556 
557 /* spBKPfactor -- sparse Bunch-Kaufman-Parlett factorisation of A in-situ
558  -- A is factored into the form P'AP = MDM' where
559  P is a permutation matrix, M lower triangular and D is block
560  diagonal with blocks of size 1 or 2
561  -- P is stored in pivot; blocks[i]==i iff D[i][i] is a block */
562 SPMAT *spBKPfactor(A,pivot,blocks,tol)
563 SPMAT *A;
564 PERM *pivot, *blocks;
565 double tol;
566 {
567  int i, j, k, l, n, onebyone=0, r;
568  int idx, idx1, idx_piv;
569  int row_num;
570  int best_deg=0, best_j, best_l, best_cost, mark_cost, deg, deg_j,
571  deg_l, ignore_deg;
572  int list_idx, list_idx2, old_list_idx;
573  SPROW *row, *r_piv, *r1_piv;
574  row_elt *e, *e1;
575  Real aii, aip1, aip1i;
576  Real det, max_j, max_l, s, t;
577  static IVEC *scan_row = IVNULL, *scan_idx = IVNULL, *col_list = IVNULL,
578  *tmp_iv = IVNULL;
579  static IVEC *deg_list = IVNULL;
580  static IVEC *orig_idx = IVNULL, *orig1_idx = IVNULL;
581  static PERM *order = PNULL;
582 
583  if ( ! A || ! pivot || ! blocks )
584  error(E_NULL,"spBKPfactor");
585  if ( A->m != A->n )
586  error(E_SQUARE,"spBKPfactor");
587  if ( A->m != pivot->size || pivot->size != blocks->size )
588  error(E_SIZES,"spBKPfactor");
589  if ( tol <= 0.0 || tol > 1.0 )
590  error(E_RANGE,"spBKPfactor");
591 
592  n = A->n;
593 
594  px_ident(pivot); px_ident(blocks);
596  ignore_deg = FALSE;
597 
598  deg_list = iv_resize(deg_list,n);
599  order = px_resize(order,n);
600  MEM_STAT_REG(deg_list,TYPE_IVEC);
602 
606  orig_idx = iv_resize(orig_idx,5);
607  orig_idx = iv_resize(orig1_idx,5);
608  orig_idx = iv_resize(tmp_iv,5);
612  MEM_STAT_REG(orig_idx,TYPE_IVEC);
613  MEM_STAT_REG(orig1_idx,TYPE_IVEC);
614  MEM_STAT_REG(tmp_iv,TYPE_IVEC);
615 
616  for ( i = 0; i < n-1; i = onebyone ? i+1 : i+2 )
617  {
618  /* now we want to use a Markowitz-style selection rule for
619  determining which rows to swap and whether to use
620  1x1 or 2x2 pivoting */
621 
622  /* get list of degrees of nodes */
623  deg_list = iv_resize(deg_list,n-i);
624  if ( ! ignore_deg )
625  for ( j = i; j < n; j++ )
626  deg_list->ive[j-i] = 0;
627  else
628  {
629  for ( j = i; j < n; j++ )
630  deg_list->ive[j-i] = 1;
631  if ( i < n )
632  deg_list->ive[0] = 0;
633  }
634  order = px_resize(order,n-i);
635  px_ident(order);
636 
637  if ( ! ignore_deg )
638  {
639  for ( j = i; j < n; j++ )
640  {
641  /* idx = sprow_idx(&(A->row[j]),j+1); */
642  /* idx = fixindex(idx); */
643  idx = 0;
644  row = &(A->row[j]);
645  e = &(row->elt[idx]);
646  /* deg_list->ive[j-i] += row->len - idx; */
647  for ( ; idx < row->len; idx++, e++ )
648  if ( e->col >= i )
649  deg_list->ive[e->col - i]++;
650  }
651  /* now deg_list[k] == degree of node k+i */
652 
653  /* now sort them into increasing order */
654  iv_sort(deg_list,order);
655  /* now deg_list[idx] == degree of node i+order[idx] */
656  }
657 
658  /* now we can chase through the nodes in order of increasing
659  degree, picking out the ones that satisfy our stability
660  criterion */
661  list_idx = 0; r = -1;
662  best_j = best_l = -1;
663  for ( deg = 0; deg <= n; deg++ )
664  {
665  Real ajj, all, ajl;
666 
667  if ( list_idx >= deg_list->dim )
668  break; /* That's all folks! */
669  old_list_idx = list_idx;
670  while ( list_idx < deg_list->dim &&
671  deg_list->ive[list_idx] <= deg )
672  {
673  j = i+order->pe[list_idx];
674  if ( j < i )
675  continue;
676  /* can we use row/col j for a 1 x 1 pivot? */
677  /* find max_j = max_{k>=i} {|A[k][j]|,|A[j][k]|} */
678  ajj = fabs(unord_get_val(A,j,j));
679  if ( ajj == 0.0 )
680  {
681  list_idx++;
682  continue; /* can't use this for 1 x 1 pivot */
683  }
684 
685  max_j = max_row_col(A,i,j,-1);
686  if ( ajj >= tol/* *alpha */ *max_j )
687  {
688  onebyone = TRUE;
689  best_j = j;
690  best_deg = deg_list->ive[list_idx];
691  break;
692  }
693  list_idx++;
694  }
695  if ( best_j >= 0 )
696  break;
697  best_cost = 2*n; /* > any possible Markowitz cost (bound) */
698  best_j = best_l = -1;
699  list_idx = old_list_idx;
700  while ( list_idx < deg_list->dim &&
701  deg_list->ive[list_idx] <= deg )
702  {
703  j = i+order->pe[list_idx];
704  ajj = fabs(unord_get_val(A,j,j));
705  for ( list_idx2 = 0; list_idx2 < list_idx; list_idx2++ )
706  {
707  deg_j = deg;
708  deg_l = deg_list->ive[list_idx2];
709  l = i+order->pe[list_idx2];
710  if ( l < i )
711  continue;
712  /* try using rows/cols (j,l) for a 2 x 2 pivot block */
713  all = fabs(unord_get_val(A,l,l));
714  ajl = ( j > l ) ? fabs(unord_get_val(A,l,j)) :
715  fabs(unord_get_val(A,j,l));
716  det = fabs(ajj*all - ajl*ajl);
717  if ( det == 0.0 )
718  continue;
719  max_j = max_row_col(A,i,j,l);
720  max_l = max_row_col(A,i,l,j);
721  if ( tol*(all*max_j+ajl*max_l) < det &&
722  tol*(ajl*max_j+ajj*max_l) < det )
723  {
724  /* acceptably stable 2 x 2 pivot */
725  /* this is actually an overestimate of the
726  Markowitz cost for choosing (j,l) */
727  mark_cost = (ajj == 0.0) ?
728  ((all == 0.0) ? deg_j+deg_l : deg_j+2*deg_l) :
729  ((all == 0.0) ? 2*deg_j+deg_l :
730  2*(deg_j+deg_l));
731  if ( mark_cost < best_cost )
732  {
733  onebyone = FALSE;
734  best_cost = mark_cost;
735  best_j = j;
736  best_l = l;
737  best_deg = deg_j;
738  }
739  }
740  }
741  list_idx++;
742  }
743  if ( best_j >= 0 )
744  break;
745  }
746 
747  if ( best_deg > (int)floor(0.8*(n-i)) )
748  ignore_deg = TRUE;
749 
750  /* now do actual interchanges */
751  if ( best_j >= 0 && onebyone )
752  {
753  bkp_interchange(A,i,best_j);
754  px_transp(pivot,i,best_j);
755  }
756  else if ( best_j >= 0 && best_l >= 0 && ! onebyone )
757  {
758  if ( best_j == i || best_j == i+1 )
759  {
760  if ( best_l == i || best_l == i+1 )
761  {
762  /* no pivoting, but must update blocks permutation */
763  px_transp(blocks,i,i+1);
764  goto dopivot;
765  }
766  bkp_interchange(A,(best_j == i) ? i+1 : i,best_l);
767  px_transp(pivot,(best_j == i) ? i+1 : i,best_l);
768  }
769  else if ( best_l == i || best_l == i+1 )
770  {
771  bkp_interchange(A,(best_l == i) ? i+1 : i,best_j);
772  px_transp(pivot,(best_l == i) ? i+1 : i,best_j);
773  }
774  else /* best_j & best_l outside i, i+1 */
775  {
776  if ( i != best_j )
777  {
778  bkp_interchange(A,i,best_j);
779  px_transp(pivot,i,best_j);
780  }
781  if ( i+1 != best_l )
782  {
783  bkp_interchange(A,i+1,best_l);
784  px_transp(pivot,i+1,best_l);
785  }
786  }
787  }
788  else /* can't pivot &/or nothing to pivot */
789  continue;
790 
791  /* update blocks permutation */
792  if ( ! onebyone )
793  px_transp(blocks,i,i+1);
794 
795  dopivot:
796  if ( onebyone )
797  {
798  int idx_j, idx_k, s_idx, s_idx2;
799  row_elt *e_ij, *e_ik;
800 
801  r_piv = &(A->row[i]);
802  idx_piv = unord_get_idx(r_piv,i);
803  /* if idx_piv < 0 then aii == 0 and no pivoting can be done;
804  -- this means that we should continue to the next iteration */
805  if ( idx_piv < 0 )
806  continue;
807  aii = r_piv->elt[idx_piv].val;
808  if ( aii == 0.0 )
809  continue;
810 
811  /* for ( j = i+1; j < n; j++ ) { ... pivot step ... } */
812  /* initialise scan_... etc for the 1 x 1 pivot */
813  scan_row = iv_resize(scan_row,r_piv->len);
814  scan_idx = iv_resize(scan_idx,r_piv->len);
815  col_list = iv_resize(col_list,r_piv->len);
816  orig_idx = iv_resize(orig_idx,r_piv->len);
817  row_num = i; s_idx = idx = 0;
818  e = &(r_piv->elt[idx]);
819  for ( idx = 0; idx < r_piv->len; idx++, e++ )
820  {
821  if ( e->col < i )
822  continue;
823  scan_row->ive[s_idx] = i;
824  scan_idx->ive[s_idx] = idx;
825  orig_idx->ive[s_idx] = idx;
826  col_list->ive[s_idx] = e->col;
827  s_idx++;
828  }
829  scan_row = iv_resize(scan_row,s_idx);
830  scan_idx = iv_resize(scan_idx,s_idx);
831  col_list = iv_resize(col_list,s_idx);
832  orig_idx = iv_resize(orig_idx,s_idx);
833 
834  order = px_resize(order,scan_row->dim);
835  px_ident(order);
837 
838  tmp_iv = iv_resize(tmp_iv,scan_row->dim);
839  for ( idx = 0; idx < order->size; idx++ )
840  tmp_iv->ive[idx] = scan_idx->ive[order->pe[idx]];
841  iv_copy(tmp_iv,scan_idx);
842  for ( idx = 0; idx < order->size; idx++ )
843  tmp_iv->ive[idx] = scan_row->ive[order->pe[idx]];
844  iv_copy(tmp_iv,scan_row);
845  for ( idx = 0; idx < scan_row->dim; idx++ )
846  tmp_iv->ive[idx] = orig_idx->ive[order->pe[idx]];
847  iv_copy(tmp_iv,orig_idx);
848 
849  /* now do actual pivot */
850  /* for ( j = i+1; j < n-1; j++ ) .... */
851 
852  for ( s_idx = 0; s_idx < scan_row->dim; s_idx++ )
853  {
854  idx_j = orig_idx->ive[s_idx];
855  if ( idx_j < 0 )
856  error(E_INTERN,"spBKPfactor");
857  e_ij = &(r_piv->elt[idx_j]);
858  j = e_ij->col;
859  if ( j < i+1 )
860  continue;
862 
863  /* compute multiplier */
864  t = e_ij->val / aii;
865 
866  /* for ( k = j; k < n; k++ ) { .... update A[j][k] .... } */
867  /* this is the row in which pivoting is done */
868  row = &(A->row[j]);
869  for ( s_idx2 = s_idx; s_idx2 < scan_row->dim; s_idx2++ )
870  {
871  idx_k = orig_idx->ive[s_idx2];
872  e_ik = &(r_piv->elt[idx_k]);
873  k = e_ik->col;
874  /* k >= j since col_list has been sorted */
875 
876  if ( scan_row->ive[s_idx2] == j )
877  { /* no fill-in -- can be done directly */
878  idx = scan_idx->ive[s_idx2];
879  /* idx = sprow_idx2(row,k,idx); */
880  row->elt[idx].val -= t*e_ik->val;
881  }
882  else
883  { /* fill-in -- insert entry & patch column */
884  int old_row, old_idx;
885  row_elt *old_e, *new_e;
886 
887  old_row = scan_row->ive[s_idx2];
888  old_idx = scan_idx->ive[s_idx2];
889  /* old_idx = sprow_idx2(&(A->row[old_row]),k,old_idx); */
890 
891  if ( old_idx < 0 )
892  error(E_INTERN,"spBKPfactor");
893  /* idx = sprow_idx(row,k); */
894  /* idx = fixindex(idx); */
895  idx = row->len;
896 
897  /* sprow_set_val(row,k,-t*e_ik->val); */
898  if ( row->len >= row->maxlen )
899  { tracecatch(sprow_xpd(row,2*row->maxlen+1,TYPE_SPMAT),
900  "spBKPfactor"); }
901 
902  row->len = idx+1;
903 
904  new_e = &(row->elt[idx]);
905  new_e->val = -t*e_ik->val;
906  new_e->col = k;
907 
908  old_e = &(A->row[old_row].elt[old_idx]);
909  new_e->nxt_row = old_e->nxt_row;
910  new_e->nxt_idx = old_e->nxt_idx;
911  old_e->nxt_row = j;
912  old_e->nxt_idx = idx;
913  }
914  }
915  e_ij->val = t;
916  }
917  }
918  else /* onebyone == FALSE */
919  { /* do 2 x 2 pivot */
920  int idx_k, idx1_k, s_idx, s_idx2;
921  int old_col;
922  row_elt *e_tmp;
923 
924  r_piv = &(A->row[i]);
925  idx_piv = unord_get_idx(r_piv,i);
926  aii = aip1i = 0.0;
927  e_tmp = r_piv->elt;
928  for ( idx_piv = 0; idx_piv < r_piv->len; idx_piv++, e_tmp++ )
929  if ( e_tmp->col == i )
930  aii = e_tmp->val;
931  else if ( e_tmp->col == i+1 )
932  aip1i = e_tmp->val;
933 
934  r1_piv = &(A->row[i+1]);
935  e_tmp = r1_piv->elt;
936  aip1 = unord_get_val(A,i+1,i+1);
937  det = aii*aip1 - aip1i*aip1i; /* Must have det < 0 */
938  if ( aii == 0.0 && aip1i == 0.0 )
939  {
940  /* error(E_RANGE,"spBKPfactor"); */
941  onebyone = TRUE;
942  continue; /* cannot pivot */
943  }
944 
945  if ( det == 0.0 )
946  {
947  if ( aii != 0.0 )
948  error(E_RANGE,"spBKPfactor");
949  onebyone = TRUE;
950  continue; /* cannot pivot */
951  }
952  aip1i = aip1i/det;
953  aii = aii/det;
954  aip1 = aip1/det;
955 
956  /* initialise scan_... etc for the 2 x 2 pivot */
957  s_idx = r_piv->len + r1_piv->len;
958  scan_row = iv_resize(scan_row,s_idx);
959  scan_idx = iv_resize(scan_idx,s_idx);
960  col_list = iv_resize(col_list,s_idx);
961  orig_idx = iv_resize(orig_idx,s_idx);
962  orig1_idx = iv_resize(orig1_idx,s_idx);
963 
964  e = r_piv->elt;
965  for ( idx = 0; idx < r_piv->len; idx++, e++ )
966  {
967  scan_row->ive[idx] = i;
968  scan_idx->ive[idx] = idx;
969  col_list->ive[idx] = e->col;
970  orig_idx->ive[idx] = idx;
971  orig1_idx->ive[idx] = -1;
972  }
973  e = r_piv->elt;
974  e1 = r1_piv->elt;
975  for ( idx = 0; idx < r1_piv->len; idx++, e1++ )
976  {
977  scan_row->ive[idx+r_piv->len] = i+1;
978  scan_idx->ive[idx+r_piv->len] = idx;
979  col_list->ive[idx+r_piv->len] = e1->col;
980  orig_idx->ive[idx+r_piv->len] = -1;
981  orig1_idx->ive[idx+r_piv->len] = idx;
982  }
983 
984  e1 = r1_piv->elt;
985  order = px_resize(order,scan_row->dim);
986  px_ident(order);
988  tmp_iv = iv_resize(tmp_iv,scan_row->dim);
989  for ( idx = 0; idx < order->size; idx++ )
990  tmp_iv->ive[idx] = scan_idx->ive[order->pe[idx]];
991  iv_copy(tmp_iv,scan_idx);
992  for ( idx = 0; idx < order->size; idx++ )
993  tmp_iv->ive[idx] = scan_row->ive[order->pe[idx]];
994  iv_copy(tmp_iv,scan_row);
995  for ( idx = 0; idx < scan_row->dim; idx++ )
996  tmp_iv->ive[idx] = orig_idx->ive[order->pe[idx]];
997  iv_copy(tmp_iv,orig_idx);
998  for ( idx = 0; idx < scan_row->dim; idx++ )
999  tmp_iv->ive[idx] = orig1_idx->ive[order->pe[idx]];
1000  iv_copy(tmp_iv,orig1_idx);
1001 
1002  s_idx = 0;
1003  old_col = -1;
1004  for ( idx = 0; idx < scan_row->dim; idx++ )
1005  {
1006  if ( col_list->ive[idx] == old_col )
1007  {
1008  if ( scan_row->ive[idx] == i )
1009  {
1010  scan_row->ive[s_idx-1] = scan_row->ive[idx];
1011  scan_idx->ive[s_idx-1] = scan_idx->ive[idx];
1012  col_list->ive[s_idx-1] = col_list->ive[idx];
1013  orig_idx->ive[s_idx-1] = orig_idx->ive[idx];
1014  orig1_idx->ive[s_idx-1] = orig1_idx->ive[idx-1];
1015  }
1016  else if ( idx > 0 )
1017  {
1018  scan_row->ive[s_idx-1] = scan_row->ive[idx-1];
1019  scan_idx->ive[s_idx-1] = scan_idx->ive[idx-1];
1020  col_list->ive[s_idx-1] = col_list->ive[idx-1];
1021  orig_idx->ive[s_idx-1] = orig_idx->ive[idx-1];
1022  orig1_idx->ive[s_idx-1] = orig1_idx->ive[idx];
1023  }
1024  }
1025  else
1026  {
1027  scan_row->ive[s_idx] = scan_row->ive[idx];
1028  scan_idx->ive[s_idx] = scan_idx->ive[idx];
1029  col_list->ive[s_idx] = col_list->ive[idx];
1030  orig_idx->ive[s_idx] = orig_idx->ive[idx];
1031  orig1_idx->ive[s_idx] = orig1_idx->ive[idx];
1032  s_idx++;
1033  }
1034  old_col = col_list->ive[idx];
1035  }
1036  scan_row = iv_resize(scan_row,s_idx);
1037  scan_idx = iv_resize(scan_idx,s_idx);
1038  col_list = iv_resize(col_list,s_idx);
1039  orig_idx = iv_resize(orig_idx,s_idx);
1040  orig1_idx = iv_resize(orig1_idx,s_idx);
1041 
1042  /* for ( j = i+2; j < n; j++ ) { .... row operation .... } */
1043  for ( s_idx = 0; s_idx < scan_row->dim; s_idx++ )
1044  {
1045  int idx_piv, idx1_piv;
1046  Real aip1j, aij, aik, aip1k;
1047  row_elt *e_ik, *e_ip1k;
1048 
1049  j = col_list->ive[s_idx];
1050  if ( j < i+2 )
1051  continue;
1053  "spBKPfactor");
1054 
1055  idx_piv = orig_idx->ive[s_idx];
1056  aij = ( idx_piv < 0 ) ? 0.0 : r_piv->elt[idx_piv].val;
1057  /* aij = ( s_idx < r_piv->len ) ? r_piv->elt[s_idx].val :
1058  0.0; */
1059  /* aij = sp_get_val(A,i,j); */
1060  idx1_piv = orig1_idx->ive[s_idx];
1061  aip1j = ( idx1_piv < 0 ) ? 0.0 : r1_piv->elt[idx1_piv].val;
1062  /* aip1j = ( s_idx < r_piv->len ) ? 0.0 :
1063  r1_piv->elt[s_idx-r_piv->len].val; */
1064  /* aip1j = sp_get_val(A,i+1,j); */
1065  s = - aip1i*aip1j + aip1*aij;
1066  t = - aip1i*aij + aii*aip1j;
1067 
1068  /* for ( k = j; k < n; k++ ) { .... update entry .... } */
1069  row = &(A->row[j]);
1070  /* set idx_k and idx1_k indices */
1071  s_idx2 = s_idx;
1072  k = col_list->ive[s_idx2];
1073  idx_k = orig_idx->ive[s_idx2];
1074  idx1_k = orig1_idx->ive[s_idx2];
1075 
1076  while ( s_idx2 < scan_row->dim )
1077  {
1078  k = col_list->ive[s_idx2];
1079  idx_k = orig_idx->ive[s_idx2];
1080  idx1_k = orig1_idx->ive[s_idx2];
1081  e_ik = ( idx_k < 0 ) ? (row_elt *)NULL :
1082  &(r_piv->elt[idx_k]);
1083  e_ip1k = ( idx1_k < 0 ) ? (row_elt *)NULL :
1084  &(r1_piv->elt[idx1_k]);
1085  aik = ( idx_k >= 0 ) ? e_ik->val : 0.0;
1086  aip1k = ( idx1_k >= 0 ) ? e_ip1k->val : 0.0;
1087  if ( scan_row->ive[s_idx2] == j )
1088  { /* no fill-in */
1089  row = &(A->row[j]);
1090  /* idx = sprow_idx(row,k); */
1091  idx = scan_idx->ive[s_idx2];
1092  if ( idx < 0 )
1093  error(E_INTERN,"spBKPfactor");
1094  row->elt[idx].val -= s*aik + t*aip1k;
1095  }
1096  else
1097  { /* fill-in -- insert entry & patch column */
1098  Real tmp;
1099  int old_row, old_idx;
1100  row_elt *old_e, *new_e;
1101 
1102  tmp = - s*aik - t*aip1k;
1103  if ( tmp != 0.0 )
1104  {
1105  row = &(A->row[j]);
1106  old_row = scan_row->ive[s_idx2];
1107  old_idx = scan_idx->ive[s_idx2];
1108 
1109  idx = row->len;
1110  if ( row->len >= row->maxlen )
1111  { tracecatch(sprow_xpd(row,2*row->maxlen+1,
1112  TYPE_SPMAT),
1113  "spBKPfactor"); }
1114 
1115  row->len = idx + 1;
1116  /* idx = sprow_idx(row,k); */
1117  new_e = &(row->elt[idx]);
1118  new_e->val = tmp;
1119  new_e->col = k;
1120 
1121  if ( old_row < 0 )
1122  error(E_INTERN,"spBKPfactor");
1123  /* old_idx = sprow_idx2(&(A->row[old_row]),
1124  k,old_idx); */
1125  old_e = &(A->row[old_row].elt[old_idx]);
1126  new_e->nxt_row = old_e->nxt_row;
1127  new_e->nxt_idx = old_e->nxt_idx;
1128  old_e->nxt_row = j;
1129  old_e->nxt_idx = idx;
1130  }
1131  }
1132 
1133  /* update idx_k, idx1_k, s_idx2 etc */
1134  s_idx2++;
1135  }
1136 
1137  /* store multipliers -- may involve fill-in (!) */
1138  /* idx = sprow_idx(r_piv,j); */
1139  idx = orig_idx->ive[s_idx];
1140  if ( idx >= 0 )
1141  {
1142  r_piv->elt[idx].val = s;
1143  }
1144  else if ( s != 0.0 )
1145  {
1146  int old_row, old_idx;
1147  row_elt *new_e, *old_e;
1148 
1149  old_row = -1; old_idx = j;
1150 
1151  if ( i > 0 )
1152  {
1153  tracecatch(chase_col(A,j,&old_row,&old_idx,i-1),
1154  "spBKPfactor");
1155  }
1156  /* sprow_set_val(r_piv,j,s); */
1157  idx = r_piv->len;
1158  if ( r_piv->len >= r_piv->maxlen )
1159  { tracecatch(sprow_xpd(r_piv,2*r_piv->maxlen+1,
1160  TYPE_SPMAT),
1161  "spBKPfactor"); }
1162 
1163  r_piv->len = idx + 1;
1164  /* idx = sprow_idx(r_piv,j); */
1165  /* if ( idx < 0 )
1166  error(E_INTERN,"spBKPfactor"); */
1167  new_e = &(r_piv->elt[idx]);
1168  new_e->val = s;
1169  new_e->col = j;
1170  if ( old_row < 0 )
1171  {
1172  new_e->nxt_row = A->start_row[j];
1173  new_e->nxt_idx = A->start_idx[j];
1174  A->start_row[j] = i;
1175  A->start_idx[j] = idx;
1176  }
1177  else
1178  {
1179  /* old_idx = sprow_idx2(&(A->row[old_row]),j,old_idx);*/
1180  if ( old_idx < 0 )
1181  error(E_INTERN,"spBKPfactor");
1182  old_e = &(A->row[old_row].elt[old_idx]);
1183  new_e->nxt_row = old_e->nxt_row;
1184  new_e->nxt_idx = old_e->nxt_idx;
1185  old_e->nxt_row = i;
1186  old_e->nxt_idx = idx;
1187  }
1188  }
1189  /* idx1 = sprow_idx(r1_piv,j); */
1190  idx1 = orig1_idx->ive[s_idx];
1191  if ( idx1 >= 0 )
1192  {
1193  r1_piv->elt[idx1].val = t;
1194  }
1195  else if ( t != 0.0 )
1196  {
1197  int old_row, old_idx;
1198  row_elt *new_e, *old_e;
1199 
1200  old_row = -1; old_idx = j;
1201  tracecatch(chase_col(A,j,&old_row,&old_idx,i),
1202  "spBKPfactor");
1203  /* sprow_set_val(r1_piv,j,t); */
1204  idx1 = r1_piv->len;
1205  if ( r1_piv->len >= r1_piv->maxlen )
1206  { tracecatch(sprow_xpd(r1_piv,2*r1_piv->maxlen+1,
1207  TYPE_SPMAT),
1208  "spBKPfactor"); }
1209 
1210  r1_piv->len = idx1 + 1;
1211  /* idx1 = sprow_idx(r1_piv,j); */
1212  /* if ( idx < 0 )
1213  error(E_INTERN,"spBKPfactor"); */
1214  new_e = &(r1_piv->elt[idx1]);
1215  new_e->val = t;
1216  new_e->col = j;
1217  if ( idx1 < 0 )
1218  error(E_INTERN,"spBKPfactor");
1219  new_e = &(r1_piv->elt[idx1]);
1220  if ( old_row < 0 )
1221  {
1222  new_e->nxt_row = A->start_row[j];
1223  new_e->nxt_idx = A->start_idx[j];
1224  A->start_row[j] = i+1;
1225  A->start_idx[j] = idx1;
1226  }
1227  else
1228  {
1229  old_idx = sprow_idx2(&(A->row[old_row]),j,old_idx);
1230  if ( old_idx < 0 )
1231  error(E_INTERN,"spBKPfactor");
1232  old_e = &(A->row[old_row].elt[old_idx]);
1233  new_e->nxt_row = old_e->nxt_row;
1234  new_e->nxt_idx = old_e->nxt_idx;
1235  old_e->nxt_row = i+1;
1236  old_e->nxt_idx = idx1;
1237  }
1238  }
1239  }
1240  }
1241  }
1242 
1243  /* now sort the rows arrays */
1244  for ( i = 0; i < A->m; i++ )
1245  qsort(A->row[i].elt,A->row[i].len,sizeof(row_elt),(int(*)())col_cmp);
1246  A->flag_col = A->flag_diag = FALSE;
1247 
1248  return A;
1249 }
1250 
1251 /* spBKPsolve -- solves A.x = b where A has been factored a la BKPfactor()
1252  -- returns x, which is created if NULL */
1253 VEC *spBKPsolve(A,pivot,block,b,x)
1254 SPMAT *A;
1255 PERM *pivot, *block;
1256 VEC *b, *x;
1257 {
1258  static VEC *tmp=VNULL; /* dummy storage needed */
1259  int i /* , j */, n, onebyone;
1260  int row_num, idx;
1261  Real a11, a12, a22, b1, b2, det, sum, *tmp_ve, tmp_diag;
1262  SPROW *r;
1263  row_elt *e;
1264 
1265  if ( ! A || ! pivot || ! block || ! b )
1266  error(E_NULL,"spBKPsolve");
1267  if ( A->m != A->n )
1268  error(E_SQUARE,"spBKPsolve");
1269  n = A->n;
1270  if ( b->dim != n || pivot->size != n || block->size != n )
1271  error(E_SIZES,"spBKPsolve");
1272  x = v_resize(x,n);
1273  tmp = v_resize(tmp,n);
1274  MEM_STAT_REG(tmp,TYPE_VEC);
1275 
1276  tmp_ve = tmp->ve;
1277 
1278  if ( ! A->flag_col )
1279  sp_col_access(A);
1280 
1281  px_vec(pivot,b,tmp);
1282  /* printf("# BKPsolve: effect of pivot: tmp =\n"); v_output(tmp); */
1283 
1284  /* solve for lower triangular part */
1285  for ( i = 0; i < n; i++ )
1286  {
1287  sum = tmp_ve[i];
1288  if ( block->pe[i] < i )
1289  {
1290  /* for ( j = 0; j < i-1; j++ )
1291  sum -= A_me[j][i]*tmp_ve[j]; */
1292  row_num = -1; idx = i;
1293  e = bump_col(A,i,&row_num,&idx);
1294  while ( row_num >= 0 && row_num < i-1 )
1295  {
1296  sum -= e->val*tmp_ve[row_num];
1297  e = bump_col(A,i,&row_num,&idx);
1298  }
1299  }
1300  else
1301  {
1302  /* for ( j = 0; j < i; j++ )
1303  sum -= A_me[j][i]*tmp_ve[j]; */
1304  row_num = -1; idx = i;
1305  e = bump_col(A,i,&row_num,&idx);
1306  while ( row_num >= 0 && row_num < i )
1307  {
1308  sum -= e->val*tmp_ve[row_num];
1309  e = bump_col(A,i,&row_num,&idx);
1310  }
1311  }
1312  tmp_ve[i] = sum;
1313  }
1314 
1315  /* printf("# BKPsolve: solving L part: tmp =\n"); v_output(tmp); */
1316  /* solve for diagonal part */
1317  for ( i = 0; i < n; i = onebyone ? i+1 : i+2 )
1318  {
1319  onebyone = ( block->pe[i] == i );
1320  if ( onebyone )
1321  {
1322  /* tmp_ve[i] /= A_me[i][i]; */
1323  tmp_diag = sp_get_val(A,i,i);
1324  if ( tmp_diag == 0.0 )
1325  error(E_SING,"spBKPsolve");
1326  tmp_ve[i] /= tmp_diag;
1327  }
1328  else
1329  {
1330  a11 = sp_get_val(A,i,i);
1331  a22 = sp_get_val(A,i+1,i+1);
1332  a12 = sp_get_val(A,i,i+1);
1333  b1 = tmp_ve[i];
1334  b2 = tmp_ve[i+1];
1335  det = a11*a22-a12*a12; /* < 0 : see BKPfactor() */
1336  if ( det == 0.0 )
1337  error(E_SING,"BKPsolve");
1338  det = 1/det;
1339  tmp_ve[i] = det*(a22*b1-a12*b2);
1340  tmp_ve[i+1] = det*(a11*b2-a12*b1);
1341  }
1342  }
1343 
1344  /* printf("# BKPsolve: solving D part: tmp =\n"); v_output(tmp); */
1345  /* solve for transpose of lower triangular part */
1346  for ( i = n-2; i >= 0; i-- )
1347  {
1348  sum = tmp_ve[i];
1349  if ( block->pe[i] > i )
1350  {
1351  /* onebyone is false */
1352  /* for ( j = i+2; j < n; j++ )
1353  sum -= A_me[i][j]*tmp_ve[j]; */
1354  if ( i+2 >= n )
1355  continue;
1356  r = &(A->row[i]);
1357  idx = sprow_idx(r,i+2);
1358  idx = fixindex(idx);
1359  e = &(r->elt[idx]);
1360  for ( ; idx < r->len; idx++, e++ )
1361  sum -= e->val*tmp_ve[e->col];
1362  }
1363  else /* onebyone */
1364  {
1365  /* for ( j = i+1; j < n; j++ )
1366  sum -= A_me[i][j]*tmp_ve[j]; */
1367  r = &(A->row[i]);
1368  idx = sprow_idx(r,i+1);
1369  idx = fixindex(idx);
1370  e = &(r->elt[idx]);
1371  for ( ; idx < r->len; idx++, e++ )
1372  sum -= e->val*tmp_ve[e->col];
1373  }
1374  tmp_ve[i] = sum;
1375  }
1376 
1377  /* printf("# BKPsolve: solving L^T part: tmp =\n");v_output(tmp); */
1378  /* and do final permutation */
1379  x = pxinv_vec(pivot,tmp,x);
1380 
1381  return x;
1382 }
1383 
1384 
1385 
short index
Definition: cabvars.h:10
static double order(void *v)
Definition: cvodeobj.cpp:239
double t
Definition: cvodeobj.cpp:59
#define TRUE
Definition: err.c:57
#define FALSE
Definition: err.c:56
#define E_SQUARE
Definition: err.h:103
#define error(err_num, fn_name)
Definition: err.h:73
#define tracecatch(ok_part, function)
Definition: err.h:166
#define E_RANGE
Definition: err.h:104
#define E_NULL
Definition: err.h:102
#define E_INTERN
Definition: err.h:111
#define E_SING
Definition: err.h:98
#define E_BOUNDS
Definition: err.h:96
#define E_SIZES
Definition: err.h:95
PERM * px_ident(PERM *px)
Definition: init.c:108
IVEC * iv_resize(IVEC *iv, int new_dim)
Definition: ivecop.c:96
IVEC * iv_sort(IVEC *x, PERM *order)
Definition: ivecop.c:359
IVEC * iv_copy(IVEC *in, IVEC *out)
Definition: ivecop.c:132
static Object ** v_resize(void *v)
Definition: ivocvect.cpp:1383
static List * done_list
Definition: kinetic.cpp:59
#define Real
Definition: machine.h:189
PERM * px_transp(PERM *px, u_int i, u_int j)
Definition: pxop.c:208
#define VNULL
Definition: matrix.h:631
#define PNULL
Definition: matrix.h:633
VEC * px_vec(PERM *, VEC *, VEC *)
PERM * px_resize(PERM *, int)
Definition: memory.c:399
#define IVNULL
Definition: matrix.h:634
VEC * pxinv_vec(PERM *, VEC *, VEC *)
#define i
Definition: md1redef.h:12
#define TYPE_PERM
Definition: meminfo.h:51
#define TYPE_VEC
Definition: meminfo.h:52
#define TYPE_IVEC
Definition: meminfo.h:53
#define MEM_STAT_REG(var, type)
Definition: meminfo.h:133
floor
Definition: extdef.h:4
fabs
Definition: extdef.h:3
#define A(i)
Definition: multisplit.cpp:63
static unsigned row
Definition: nonlin.cpp:85
int const size_t const size_t n
Definition: nrngsl.h:11
for(i=0;i< n;i++)
if(status)
size_t j
static philox4x32_key_t k
Definition: nrnran123.cpp:11
#define e
Definition: passive0.cpp:22
row_elt * bump_col(SPMAT *, int, int *, int *)
Definition: spswap.c:258
void scan_to(SPMAT *, IVEC *, IVEC *, IVEC *, int)
Definition: spswap.c:45
row_elt * chase_col(SPMAT *, int, int *, int *, int)
Definition: spswap.c:134
row_elt * chase_past(SPMAT *, int, int *, int *, int)
Definition: spswap.c:210
double sp_get_val(SPMAT *A, int i, int j)
Definition: sparse.c:45
SPMAT * sp_diag_access(SPMAT *A)
Definition: sparse.c:417
SPMAT * sp_col_access(SPMAT *A)
Definition: sparse.c:376
SPROW * sprow_xpd(SPROW *r, int n, int type)
#define fixindex(idx)
Definition: sparse.h:186
int sprow_idx(SPROW *, int)
Definition: sprow.c:75
#define sprow_idx2(r, c, hint)
Definition: sparse.h:78
row_elt * bkp_bump_col(SPMAT *A, int j, int *row, int *idx)
Definition: spbkp.c:203
int chk_col_access(SPMAT *A)
Definition: spbkp.c:512
SPMAT * bkp_interchange(SPMAT *A, int i1, int i2)
Definition: spbkp.c:232
SPMAT * spBKPfactor(SPMAT *A, PERM *pivot, PERM *blocks, double tol)
Definition: spbkp.c:562
double unord_get_val(SPMAT *A, int i, int j)
Definition: spbkp.c:77
static SPMAT * bkp_swap_elt(SPMAT *A, int i1, int j1, int idx1, int i2, int j2, int idx2)
Definition: spbkp.c:100
static int nonzeros(SPMAT *A)
Definition: spbkp.c:496
int iv_min(IVEC *iv, int *index)
Definition: spbkp.c:414
static double max_row_col(SPMAT *A, int i, int j, int l)
Definition: spbkp.c:444
int unord_get_idx(SPROW *r, int j)
Definition: spbkp.c:57
static int col_cmp(row_elt *e1, row_elt *e2)
Definition: spbkp.c:551
VEC * spBKPsolve(SPMAT *A, PERM *pivot, PERM *block, VEC *b, VEC *x)
Definition: spbkp.c:1253
static char rcsid[]
Definition: spbkp.c:35
static int * col_list
Definition: spchfctr.c:147
static int * scan_row
Definition: spchfctr.c:146
static int * scan_idx
Definition: spchfctr.c:146
#define cnt
Definition: spt2queue.cpp:19
#define NULL
Definition: sptree.h:16
Definition: matrix.h:92
int * ive
Definition: matrix.h:94
u_int dim
Definition: matrix.h:93
Definition: matrix.h:87
u_int * pe
Definition: matrix.h:88
u_int size
Definition: matrix.h:88
Definition: sparse.h:54
Definition: sparse.h:49
int maxlen
Definition: sparse.h:50
int len
Definition: sparse.h:50
row_elt * elt
Definition: sparse.h:51
Definition: matrix.h:67
u_int dim
Definition: matrix.h:68
Real * ve
Definition: matrix.h:69
Definition: sparse.h:44
int nxt_row
Definition: sparse.h:45
int nxt_idx
Definition: sparse.h:45
int col
Definition: sparse.h:45
Real val
Definition: sparse.h:46