xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision 517f4e3460cd8c7e73c19b8bf533f1efe47f30a7)
1 /*
2               This file creating by running f2c
3             linpack. this version dated 08/14/78
4       cleve moler, university of new mexico, argonne national lab.
5 
6       Computes the inverse of a matrix given its factors and pivots
7     calculated by PetscLINPACKgefa(). Performed in-place for an n by n
8     dense matrix.
9 
10        Used by the sparse factorization routines in
11      src/mat/impls/baij/seq
12 
13 */
14 
15 #include <petscsys.h>
16 #include <petsc/private/kernels/blockinvert.h>
17 
18 PETSC_INTERN PetscErrorCode PetscLINPACKgedi(MatScalar *a, PetscInt n, PetscInt *ipvt, MatScalar *work)
19 {
20   PetscInt   i__2, kb, kp1, nm1, i, j, k, l, ll, kn, knp1, jn1;
21   MatScalar *aa, *ax, *ay, tmp;
22   MatScalar  t;
23 
24   PetscFunctionBegin;
25   --work;
26   --ipvt;
27   a -= n + 1;
28 
29   /*     compute inverse(u) */
30 
31   for (k = 1; k <= n; ++k) {
32     kn      = k * n;
33     knp1    = kn + k;
34     a[knp1] = 1.0 / a[knp1];
35     t       = -a[knp1];
36     i__2    = k - 1;
37     aa      = &a[1 + kn];
38     for (ll = 0; ll < i__2; ll++) aa[ll] *= t;
39     kp1 = k + 1;
40     if (n < kp1) continue;
41     ax = aa;
42     for (j = kp1; j <= n; ++j) {
43       jn1        = j * n;
44       t          = a[k + jn1];
45       a[k + jn1] = 0.;
46       ay         = &a[1 + jn1];
47       for (ll = 0; ll < k; ll++) ay[ll] += t * ax[ll];
48     }
49   }
50 
51   /*    form inverse(u)*inverse(l) */
52 
53   nm1 = n - 1;
54   if (nm1 < 1) PetscFunctionReturn(PETSC_SUCCESS);
55   for (kb = 1; kb <= nm1; ++kb) {
56     k   = n - kb;
57     kn  = k * n;
58     kp1 = k + 1;
59     aa  = a + kn;
60     for (i = kp1; i <= n; ++i) {
61       work[i] = aa[i];
62       aa[i]   = 0.;
63     }
64     for (j = kp1; j <= n; ++j) {
65       t  = work[j];
66       ax = &a[j * n + 1];
67       ay = &a[kn + 1];
68       for (ll = 0; ll < n; ll++) ay[ll] += t * ax[ll];
69     }
70     l = ipvt[k];
71     if (l != k) {
72       ax = &a[kn + 1];
73       ay = &a[l * n + 1];
74       for (ll = 0; ll < n; ll++) {
75         tmp    = ax[ll];
76         ax[ll] = ay[ll];
77         ay[ll] = tmp;
78       }
79     }
80   }
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83