Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
amesos_klu_l.c
Go to the documentation of this file.
1 /* ========================================================================== */
2 /* === klu ================================================================== */
3 /* ========================================================================== */
4 
5 /* KLU: factorizes P*A into L*U, using the Gilbert-Peierls algorithm [1], with
6  * optional symmetric pruning by Eisenstat and Liu [2]. The code is by Tim
7  * Davis. This algorithm is what appears as the default sparse LU routine in
8  * MATLAB version 6.0, and still appears in MATLAB 6.5 as [L,U,P] = lu (A).
9  * Note that no column ordering is provided (see COLAMD or AMD for suitable
10  * orderings). SuperLU is based on this algorithm, except that it adds the
11  * use of dense matrix operations on "supernodes" (adjacent columns with
12  * identical). This code doesn't use supernodes, thus its name ("Kent" LU,
13  * as in "Clark Kent", in contrast with Super-LU...). This algorithm is slower
14  * than SuperLU and UMFPACK for large matrices with lots of nonzeros in their
15  * factors (such as for most finite-element problems). However, for matrices
16  * with very sparse LU factors, this algorithm is typically faster than both
17  * SuperLU and UMFPACK, since in this case there is little chance to exploit
18  * dense matrix kernels (the BLAS).
19  *
20  * Only one block of A is factorized, in the BTF form. The input n is the
21  * size of the block; k1 is the first row and column in the block.
22  *
23  * NOTE: no error checking is done on the inputs. This version is not meant to
24  * be called directly by the user. Use klu_factor instead.
25  *
26  * No fill-reducing ordering is provided. The ordering quality of
27  * klu_kernel_factor is the responsibility of the caller. The input A must
28  * pre-permuted to reduce fill-in, or fill-reducing input permutation Q must
29  * be provided.
30  *
31  * The input matrix A must be in compressed-column form, with either sorted
32  * or unsorted row indices. Row indices for column j of A is in
33  * Ai [Ap [j] ... Ap [j+1]-1] and the same range of indices in Ax holds the
34  * numerical values. No duplicate entries are allowed.
35  *
36  * Copyright 2004-2007, Tim Davis. All rights reserved. See the README
37  * file for details on permitted use. Note that no code from The MathWorks,
38  * Inc, or from SuperLU, or from any other source appears here. The code is
39  * written from scratch, from the algorithmic description in Gilbert & Peierls'
40  * and Eisenstat & Liu's journal papers [1,2].
41  *
42  * If an input permutation Q is provided, the factorization L*U = A (P,Q)
43  * is computed, where P is determined by partial pivoting, and Q is the input
44  * ordering. If the pivot tolerance is less than 1, the "diagonal" entry that
45  * KLU attempts to choose is the diagonal of A (Q,Q). In other words, the
46  * input permutation is applied symmetrically to the input matrix. The output
47  * permutation P includes both the partial pivoting ordering and the input
48  * permutation. If Q is NULL, then it is assumed to be the identity
49  * permutation. Q is not modified.
50  *
51  * [1] Gilbert, J. R. and Peierls, T., "Sparse Partial Pivoting in Time
52  * Proportional to Arithmetic Operations," SIAM J. Sci. Stat. Comp.,
53  * vol 9, pp. 862-874, 1988.
54  * [2] Eisenstat, S. C. and Liu, J. W. H., "Exploiting Structural Symmetry in
55  * Unsymmetric Sparse Symbolic Factorization," SIAM J. Matrix Analysis &
56  * Applic., vol 13, pp. 202-211, 1992.
57  */
58 
59 /* ========================================================================== */
60 
61 /* This file should make the long int version of KLU */
62 #define DLONG 1
63 
64 #include "amesos_klu_internal.h"
65 
66 size_t KLU_kernel_factor /* 0 if failure, size of LU if OK */
67 (
68  /* inputs, not modified */
69  Int n, /* A is n-by-n. n must be > 0. */
70  Int Ap [ ], /* size n+1, column pointers for A */
71  Int Ai [ ], /* size nz = Ap [n], row indices for A */
72  Entry Ax [ ], /* size nz, values of A */
73  Int Q [ ], /* size n, optional column permutation */
74  double Lsize, /* estimate of number of nonzeros in L */
75 
76  /* outputs, not defined on input */
77  Unit **p_LU, /* row indices and values of L and U */
78  Entry Udiag [ ], /* size n, diagonal of U */
79  Int Llen [ ], /* size n, column length of L */
80  Int Ulen [ ], /* size n, column length of U */
81  Int Lip [ ], /* size n, column pointers for L */
82  Int Uip [ ], /* size n, column pointers for U */
83  Int P [ ], /* row permutation, size n */
84  Int *lnz, /* size of L */
85  Int *unz, /* size of U */
86 
87  /* workspace, undefined on input */
88  Entry *X, /* size n double's, zero on output */
89  Int *Work, /* size 5n Int's */
90 
91  /* inputs, not modified on output */
92  Int k1, /* the block of A is from k1 to k2-1 */
93  Int PSinv [ ], /* inverse of P from symbolic factorization */
94  double Rs [ ], /* scale factors for A */
95 
96  /* inputs, modified on output */
97  Int Offp [ ], /* off-diagonal matrix (modified by this routine) */
98  Int Offi [ ],
99  Entry Offx [ ],
100  /* --------------- */
101  KLU_common *Common
102 )
103 {
104  double maxlnz, dunits ;
105  Unit *LU ;
106  Int *Pinv, *Lpend, *Stack, *Flag, *Ap_pos, *W ;
107  Int lsize, usize, anz, ok ;
108  size_t lusize ;
109  ASSERT (Common != NULL) ;
110 
111  /* ---------------------------------------------------------------------- */
112  /* get control parameters, or use defaults */
113  /* ---------------------------------------------------------------------- */
114 
115  n = MAX (1, n) ;
116  anz = Ap [n+k1] - Ap [k1] ;
117 
118  if (Lsize <= 0)
119  {
120  Lsize = -Lsize ;
121  Lsize = MAX (Lsize, 1.0) ;
122  lsize = Lsize * anz + n ;
123  }
124  else
125  {
126  lsize = Lsize ;
127  }
128 
129  usize = lsize ;
130 
131  lsize = MAX (n+1, lsize) ;
132  usize = MAX (n+1, usize) ;
133 
134  maxlnz = (((double) n) * ((double) n) + ((double) n)) / 2. ;
135  maxlnz = MIN (maxlnz, ((double) INT_MAX)) ;
136  lsize = MIN (maxlnz, lsize) ;
137  usize = MIN (maxlnz, usize) ;
138 
139  PRINTF (("Welcome to klu: n %d anz %d k1 %d lsize %d usize %d maxlnz %g\n",
140  n, anz, k1, lsize, usize, maxlnz)) ;
141 
142  /* ---------------------------------------------------------------------- */
143  /* allocate workspace and outputs */
144  /* ---------------------------------------------------------------------- */
145 
146  /* return arguments are not yet assigned */
147  *p_LU = (Unit *) NULL ;
148 
149  /* these computations are safe from size_t overflow */
150  W = Work ;
151  Pinv = (Int *) W ; W += n ;
152  Stack = (Int *) W ; W += n ;
153  Flag = (Int *) W ; W += n ;
154  Lpend = (Int *) W ; W += n ;
155  Ap_pos = (Int *) W ; W += n ;
156 
157  dunits = DUNITS (Int, lsize) + DUNITS (Entry, lsize) +
158  DUNITS (Int, usize) + DUNITS (Entry, usize) ;
159  lusize = (size_t) dunits ;
160  ok = !INT_OVERFLOW (dunits) ;
161  LU = ok ? KLU_malloc (lusize, sizeof (Unit), Common) : NULL ;
162  if (LU == NULL)
163  {
164  /* out of memory, or problem too large */
165  Common->status = KLU_OUT_OF_MEMORY ;
166  lusize = 0 ;
167  return (lusize) ;
168  }
169 
170  /* ---------------------------------------------------------------------- */
171  /* factorize */
172  /* ---------------------------------------------------------------------- */
173 
174  /* with pruning, and non-recursive depth-first-search */
175  lusize = KLU_kernel (n, Ap, Ai, Ax, Q, lusize,
176  Pinv, P, &LU, Udiag, Llen, Ulen, Lip, Uip, lnz, unz,
177  X, Stack, Flag, Ap_pos, Lpend,
178  k1, PSinv, Rs, Offp, Offi, Offx, Common) ;
179 
180  /* ---------------------------------------------------------------------- */
181  /* return LU factors, or return nothing if an error occurred */
182  /* ---------------------------------------------------------------------- */
183 
184  if (Common->status < KLU_OK)
185  {
186  LU = KLU_free (LU, lusize, sizeof (Unit), Common) ;
187  lusize = 0 ;
188  }
189  *p_LU = LU ;
190  PRINTF ((" in klu noffdiag %d\n", Common->noffdiag)) ;
191  return (lusize) ;
192 }
193 
194 
195 /* ========================================================================== */
196 /* === KLU_lsolve =========================================================== */
197 /* ========================================================================== */
198 
199 /* Solve Lx=b. Assumes L is unit lower triangular and where the unit diagonal
200  * entry is NOT stored. Overwrites B with the solution X. B is n-by-nrhs
201  * and is stored in ROW form with row dimension nrhs. nrhs must be in the
202  * range 1 to 4. */
203 
204 void KLU_lsolve
205 (
206  /* inputs, not modified: */
207  Int n,
208  Int Lip [ ],
209  Int Llen [ ],
210  Unit LU [ ],
211  Int nrhs,
212  /* right-hand-side on input, solution to Lx=b on output */
213  Entry X [ ]
214 )
215 {
216  Entry x [4], lik ;
217  Int *Li ;
218  Entry *Lx ;
219  Int k, p, len, i ;
220 
221  switch (nrhs)
222  {
223 
224  case 1:
225  for (k = 0 ; k < n ; k++)
226  {
227  x [0] = X [k] ;
228  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
229  /* unit diagonal of L is not stored*/
230  for (p = 0 ; p < len ; p++)
231  {
232  /* X [Li [p]] -= Lx [p] * x [0] ; */
233  MULT_SUB (X [Li [p]], Lx [p], x [0]) ;
234  }
235  }
236  break ;
237 
238  case 2:
239 
240  for (k = 0 ; k < n ; k++)
241  {
242  x [0] = X [2*k ] ;
243  x [1] = X [2*k + 1] ;
244  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
245  for (p = 0 ; p < len ; p++)
246  {
247  i = Li [p] ;
248  lik = Lx [p] ;
249  MULT_SUB (X [2*i], lik, x [0]) ;
250  MULT_SUB (X [2*i + 1], lik, x [1]) ;
251  }
252  }
253  break ;
254 
255  case 3:
256 
257  for (k = 0 ; k < n ; k++)
258  {
259  x [0] = X [3*k ] ;
260  x [1] = X [3*k + 1] ;
261  x [2] = X [3*k + 2] ;
262  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
263  for (p = 0 ; p < len ; p++)
264  {
265  i = Li [p] ;
266  lik = Lx [p] ;
267  MULT_SUB (X [3*i], lik, x [0]) ;
268  MULT_SUB (X [3*i + 1], lik, x [1]) ;
269  MULT_SUB (X [3*i + 2], lik, x [2]) ;
270  }
271  }
272  break ;
273 
274  case 4:
275 
276  for (k = 0 ; k < n ; k++)
277  {
278  x [0] = X [4*k ] ;
279  x [1] = X [4*k + 1] ;
280  x [2] = X [4*k + 2] ;
281  x [3] = X [4*k + 3] ;
282  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
283  for (p = 0 ; p < len ; p++)
284  {
285  i = Li [p] ;
286  lik = Lx [p] ;
287  MULT_SUB (X [4*i], lik, x [0]) ;
288  MULT_SUB (X [4*i + 1], lik, x [1]) ;
289  MULT_SUB (X [4*i + 2], lik, x [2]) ;
290  MULT_SUB (X [4*i + 3], lik, x [3]) ;
291  }
292  }
293  break ;
294 
295  }
296 }
297 
298 /* ========================================================================== */
299 /* === KLU_usolve =========================================================== */
300 /* ========================================================================== */
301 
302 /* Solve Ux=b. Assumes U is non-unit upper triangular and where the diagonal
303  * entry is NOT stored. Overwrites B with the solution X. B is n-by-nrhs
304  * and is stored in ROW form with row dimension nrhs. nrhs must be in the
305  * range 1 to 4. */
306 
307 void KLU_usolve
308 (
309  /* inputs, not modified: */
310  Int n,
311  Int Uip [ ],
312  Int Ulen [ ],
313  Unit LU [ ],
314  Entry Udiag [ ],
315  Int nrhs,
316  /* right-hand-side on input, solution to Ux=b on output */
317  Entry X [ ]
318 )
319 {
320  Entry x [4], uik, ukk ;
321  Int *Ui ;
322  Entry *Ux ;
323  Int k, p, len, i ;
324 
325  switch (nrhs)
326  {
327 
328  case 1:
329 
330  for (k = n-1 ; k >= 0 ; k--)
331  {
332  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
333  /* x [0] = X [k] / Udiag [k] ; */
334  DIV (x [0], X [k], Udiag [k]) ;
335  X [k] = x [0] ;
336  for (p = 0 ; p < len ; p++)
337  {
338  /* X [Ui [p]] -= Ux [p] * x [0] ; */
339  MULT_SUB (X [Ui [p]], Ux [p], x [0]) ;
340 
341  }
342  }
343 
344  break ;
345 
346  case 2:
347 
348  for (k = n-1 ; k >= 0 ; k--)
349  {
350  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
351  ukk = Udiag [k] ;
352  /* x [0] = X [2*k ] / ukk ;
353  x [1] = X [2*k + 1] / ukk ; */
354  DIV (x [0], X [2*k], ukk) ;
355  DIV (x [1], X [2*k + 1], ukk) ;
356 
357  X [2*k ] = x [0] ;
358  X [2*k + 1] = x [1] ;
359  for (p = 0 ; p < len ; p++)
360  {
361  i = Ui [p] ;
362  uik = Ux [p] ;
363  /* X [2*i ] -= uik * x [0] ;
364  X [2*i + 1] -= uik * x [1] ; */
365  MULT_SUB (X [2*i], uik, x [0]) ;
366  MULT_SUB (X [2*i + 1], uik, x [1]) ;
367  }
368  }
369 
370  break ;
371 
372  case 3:
373 
374  for (k = n-1 ; k >= 0 ; k--)
375  {
376  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
377  ukk = Udiag [k] ;
378 
379  DIV (x [0], X [3*k], ukk) ;
380  DIV (x [1], X [3*k + 1], ukk) ;
381  DIV (x [2], X [3*k + 2], ukk) ;
382 
383  X [3*k ] = x [0] ;
384  X [3*k + 1] = x [1] ;
385  X [3*k + 2] = x [2] ;
386  for (p = 0 ; p < len ; p++)
387  {
388  i = Ui [p] ;
389  uik = Ux [p] ;
390  MULT_SUB (X [3*i], uik, x [0]) ;
391  MULT_SUB (X [3*i + 1], uik, x [1]) ;
392  MULT_SUB (X [3*i + 2], uik, x [2]) ;
393  }
394  }
395 
396  break ;
397 
398  case 4:
399 
400  for (k = n-1 ; k >= 0 ; k--)
401  {
402  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
403  ukk = Udiag [k] ;
404 
405  DIV (x [0], X [4*k], ukk) ;
406  DIV (x [1], X [4*k + 1], ukk) ;
407  DIV (x [2], X [4*k + 2], ukk) ;
408  DIV (x [3], X [4*k + 3], ukk) ;
409 
410  X [4*k ] = x [0] ;
411  X [4*k + 1] = x [1] ;
412  X [4*k + 2] = x [2] ;
413  X [4*k + 3] = x [3] ;
414  for (p = 0 ; p < len ; p++)
415  {
416  i = Ui [p] ;
417  uik = Ux [p] ;
418 
419  MULT_SUB (X [4*i], uik, x [0]) ;
420  MULT_SUB (X [4*i + 1], uik, x [1]) ;
421  MULT_SUB (X [4*i + 2], uik, x [2]) ;
422  MULT_SUB (X [4*i + 3], uik, x [3]) ;
423  }
424  }
425 
426  break ;
427 
428  }
429 }
430 
431 
432 /* ========================================================================== */
433 /* === KLU_ltsolve ========================================================== */
434 /* ========================================================================== */
435 
436 /* Solve L'x=b. Assumes L is unit lower triangular and where the unit diagonal
437  * entry is NOT stored. Overwrites B with the solution X. B is n-by-nrhs
438  * and is stored in ROW form with row dimension nrhs. nrhs must in the
439  * range 1 to 4. */
440 
441 void KLU_ltsolve
442 (
443  /* inputs, not modified: */
444  Int n,
445  Int Lip [ ],
446  Int Llen [ ],
447  Unit LU [ ],
448  Int nrhs,
449 #ifdef COMPLEX
450  Int conj_solve,
451 #endif
452  /* right-hand-side on input, solution to L'x=b on output */
453  Entry X [ ]
454 )
455 {
456  Entry x [4], lik ;
457  Int *Li ;
458  Entry *Lx ;
459  Int k, p, len, i ;
460 
461  switch (nrhs)
462  {
463 
464  case 1:
465 
466  for (k = n-1 ; k >= 0 ; k--)
467  {
468  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
469  x [0] = X [k] ;
470  for (p = 0 ; p < len ; p++)
471  {
472 #ifdef COMPLEX
473  if (conj_solve)
474  {
475  /* x [0] -= CONJ (Lx [p]) * X [Li [p]] ; */
476  MULT_SUB_CONJ (x [0], X [Li [p]], Lx [p]) ;
477  }
478  else
479 #endif
480  {
481  /*x [0] -= Lx [p] * X [Li [p]] ;*/
482  MULT_SUB (x [0], Lx [p], X [Li [p]]) ;
483  }
484  }
485  X [k] = x [0] ;
486  }
487  break ;
488 
489  case 2:
490 
491  for (k = n-1 ; k >= 0 ; k--)
492  {
493  x [0] = X [2*k ] ;
494  x [1] = X [2*k + 1] ;
495  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
496  for (p = 0 ; p < len ; p++)
497  {
498  i = Li [p] ;
499 #ifdef COMPLEX
500  if (conj_solve)
501  {
502  CONJ (lik, Lx [p]) ;
503  }
504  else
505 #endif
506  {
507  lik = Lx [p] ;
508  }
509  MULT_SUB (x [0], lik, X [2*i]) ;
510  MULT_SUB (x [1], lik, X [2*i + 1]) ;
511  }
512  X [2*k ] = x [0] ;
513  X [2*k + 1] = x [1] ;
514  }
515  break ;
516 
517  case 3:
518 
519  for (k = n-1 ; k >= 0 ; k--)
520  {
521  x [0] = X [3*k ] ;
522  x [1] = X [3*k + 1] ;
523  x [2] = X [3*k + 2] ;
524  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
525  for (p = 0 ; p < len ; p++)
526  {
527  i = Li [p] ;
528 #ifdef COMPLEX
529  if (conj_solve)
530  {
531  CONJ (lik, Lx [p]) ;
532  }
533  else
534 #endif
535  {
536  lik = Lx [p] ;
537  }
538  MULT_SUB (x [0], lik, X [3*i]) ;
539  MULT_SUB (x [1], lik, X [3*i + 1]) ;
540  MULT_SUB (x [2], lik, X [3*i + 2]) ;
541  }
542  X [3*k ] = x [0] ;
543  X [3*k + 1] = x [1] ;
544  X [3*k + 2] = x [2] ;
545  }
546  break ;
547 
548  case 4:
549 
550  for (k = n-1 ; k >= 0 ; k--)
551  {
552  x [0] = X [4*k ] ;
553  x [1] = X [4*k + 1] ;
554  x [2] = X [4*k + 2] ;
555  x [3] = X [4*k + 3] ;
556  GET_POINTER (LU, Lip, Llen, Li, Lx, k, len) ;
557  for (p = 0 ; p < len ; p++)
558  {
559  i = Li [p] ;
560 #ifdef COMPLEX
561  if (conj_solve)
562  {
563  CONJ (lik, Lx [p]) ;
564  }
565  else
566 #endif
567  {
568  lik = Lx [p] ;
569  }
570  MULT_SUB (x [0], lik, X [4*i]) ;
571  MULT_SUB (x [1], lik, X [4*i + 1]) ;
572  MULT_SUB (x [2], lik, X [4*i + 2]) ;
573  MULT_SUB (x [3], lik, X [4*i + 3]) ;
574  }
575  X [4*k ] = x [0] ;
576  X [4*k + 1] = x [1] ;
577  X [4*k + 2] = x [2] ;
578  X [4*k + 3] = x [3] ;
579  }
580  break ;
581  }
582 }
583 
584 
585 /* ========================================================================== */
586 /* === KLU_utsolve ========================================================== */
587 /* ========================================================================== */
588 
589 /* Solve U'x=b. Assumes U is non-unit upper triangular and where the diagonal
590  * entry is stored (and appears last in each column of U). Overwrites B
591  * with the solution X. B is n-by-nrhs and is stored in ROW form with row
592  * dimension nrhs. nrhs must be in the range 1 to 4. */
593 
594 void KLU_utsolve
595 (
596  /* inputs, not modified: */
597  Int n,
598  Int Uip [ ],
599  Int Ulen [ ],
600  Unit LU [ ],
601  Entry Udiag [ ],
602  Int nrhs,
603 #ifdef COMPLEX
604  Int conj_solve,
605 #endif
606  /* right-hand-side on input, solution to Ux=b on output */
607  Entry X [ ]
608 )
609 {
610  Entry x [4], uik, ukk ;
611  Int k, p, len, i ;
612  Int *Ui ;
613  Entry *Ux ;
614 
615  switch (nrhs)
616  {
617 
618  case 1:
619 
620  for (k = 0 ; k < n ; k++)
621  {
622  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
623  x [0] = X [k] ;
624  for (p = 0 ; p < len ; p++)
625  {
626 #ifdef COMPLEX
627  if (conj_solve)
628  {
629  /* x [0] -= CONJ (Ux [p]) * X [Ui [p]] ; */
630  MULT_SUB_CONJ (x [0], X [Ui [p]], Ux [p]) ;
631  }
632  else
633 #endif
634  {
635  /* x [0] -= Ux [p] * X [Ui [p]] ; */
636  MULT_SUB (x [0], Ux [p], X [Ui [p]]) ;
637  }
638  }
639 #ifdef COMPLEX
640  if (conj_solve)
641  {
642  CONJ (ukk, Udiag [k]) ;
643  }
644  else
645 #endif
646  {
647  ukk = Udiag [k] ;
648  }
649  DIV (X [k], x [0], ukk) ;
650  }
651  break ;
652 
653  case 2:
654 
655  for (k = 0 ; k < n ; k++)
656  {
657  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
658  x [0] = X [2*k ] ;
659  x [1] = X [2*k + 1] ;
660  for (p = 0 ; p < len ; p++)
661  {
662  i = Ui [p] ;
663 #ifdef COMPLEX
664  if (conj_solve)
665  {
666  CONJ (uik, Ux [p]) ;
667  }
668  else
669 #endif
670  {
671  uik = Ux [p] ;
672  }
673  MULT_SUB (x [0], uik, X [2*i]) ;
674  MULT_SUB (x [1], uik, X [2*i + 1]) ;
675  }
676 #ifdef COMPLEX
677  if (conj_solve)
678  {
679  CONJ (ukk, Udiag [k]) ;
680  }
681  else
682 #endif
683  {
684  ukk = Udiag [k] ;
685  }
686  DIV (X [2*k], x [0], ukk) ;
687  DIV (X [2*k + 1], x [1], ukk) ;
688  }
689  break ;
690 
691  case 3:
692 
693  for (k = 0 ; k < n ; k++)
694  {
695  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
696  x [0] = X [3*k ] ;
697  x [1] = X [3*k + 1] ;
698  x [2] = X [3*k + 2] ;
699  for (p = 0 ; p < len ; p++)
700  {
701  i = Ui [p] ;
702 #ifdef COMPLEX
703  if (conj_solve)
704  {
705  CONJ (uik, Ux [p]) ;
706  }
707  else
708 #endif
709  {
710  uik = Ux [p] ;
711  }
712  MULT_SUB (x [0], uik, X [3*i]) ;
713  MULT_SUB (x [1], uik, X [3*i + 1]) ;
714  MULT_SUB (x [2], uik, X [3*i + 2]) ;
715  }
716 #ifdef COMPLEX
717  if (conj_solve)
718  {
719  CONJ (ukk, Udiag [k]) ;
720  }
721  else
722 #endif
723  {
724  ukk = Udiag [k] ;
725  }
726  DIV (X [3*k], x [0], ukk) ;
727  DIV (X [3*k + 1], x [1], ukk) ;
728  DIV (X [3*k + 2], x [2], ukk) ;
729  }
730  break ;
731 
732  case 4:
733 
734  for (k = 0 ; k < n ; k++)
735  {
736  GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, len) ;
737  x [0] = X [4*k ] ;
738  x [1] = X [4*k + 1] ;
739  x [2] = X [4*k + 2] ;
740  x [3] = X [4*k + 3] ;
741  for (p = 0 ; p < len ; p++)
742  {
743  i = Ui [p] ;
744 #ifdef COMPLEX
745  if (conj_solve)
746  {
747  CONJ (uik, Ux [p]) ;
748  }
749  else
750 #endif
751  {
752  uik = Ux [p] ;
753  }
754  MULT_SUB (x [0], uik, X [4*i]) ;
755  MULT_SUB (x [1], uik, X [4*i + 1]) ;
756  MULT_SUB (x [2], uik, X [4*i + 2]) ;
757  MULT_SUB (x [3], uik, X [4*i + 3]) ;
758  }
759 #ifdef COMPLEX
760  if (conj_solve)
761  {
762  CONJ (ukk, Udiag [k]) ;
763  }
764  else
765 #endif
766  {
767  ukk = Udiag [k] ;
768  }
769  DIV (X [4*k], x [0], ukk) ;
770  DIV (X [4*k + 1], x [1], ukk) ;
771  DIV (X [4*k + 2], x [2], ukk) ;
772  DIV (X [4*k + 3], x [3], ukk) ;
773  }
774  break ;
775  }
776 }
#define GET_POINTER(LU, Xip, Xlen, Xi, Xx, k, xlen)
#define CONJ(a, x)
#define COMPLEX
#define Int
#define KLU_kernel
#define P(k)
void KLU_utsolve(Int n, Int Uip[], Int Ulen[], Unit LU[], Entry Udiag[], Int nrhs, Entry X[])
Definition: amesos_klu_l.c:595
double Unit
#define MAX(a, b)
#define KLU_OUT_OF_MEMORY
#define NULL
#define ASSERT(expression)
#define MULT_SUB(c, a, b)
#define KLU_malloc
#define MULT_SUB_CONJ(c, a, b)
#define INT_OVERFLOW(x)
size_t KLU_kernel_factor(Int n, Int Ap[], Int Ai[], Entry Ax[], Int Q[], double Lsize, Unit **p_LU, Entry Udiag[], Int Llen[], Int Ulen[], Int Lip[], Int Uip[], Int P[], Int *lnz, Int *unz, Entry *X, Int *Work, Int k1, Int PSinv[], double Rs[], Int Offp[], Int Offi[], Entry Offx[], KLU_common *Common)
Definition: amesos_klu_l.c:67
#define DIV(c, a, b)
#define KLU_common
#define DUNITS(type, n)
#define MIN(a, b)
void KLU_lsolve(Int n, Int Lip[], Int Llen[], Unit LU[], Int nrhs, Entry X[])
Definition: amesos_klu_l.c:205
void KLU_usolve(Int n, Int Uip[], Int Ulen[], Unit LU[], Entry Udiag[], Int nrhs, Entry X[])
Definition: amesos_klu_l.c:308
#define Entry
#define PRINTF(params)
void KLU_ltsolve(Int n, Int Lip[], Int Llen[], Unit LU[], Int nrhs, Entry X[])
Definition: amesos_klu_l.c:442
int n
#define KLU_free
#define KLU_OK