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