xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
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 #include <petsc/private/kernels/blockinvert.h>
18 
19 PETSC_INTERN PetscErrorCode PetscLINPACKgedi(MatScalar *a,PetscInt n,PetscInt *ipvt,MatScalar *work)
20 {
21   PetscInt  i__2,kb,kp1,nm1,i,j,k,l,ll,kn,knp1,jn1;
22   MatScalar *aa,*ax,*ay,tmp;
23   MatScalar t;
24 
25   PetscFunctionBegin;
26   --work;
27   --ipvt;
28   a -= n + 1;
29 
30   /*     compute inverse(u) */
31 
32   for (k = 1; k <= n; ++k) {
33     kn      = k*n;
34     knp1    = kn + k;
35     a[knp1] = 1.0 / a[knp1];
36     t       = -a[knp1];
37     i__2    = k - 1;
38     aa      = &a[1 + kn];
39     for (ll=0; ll<i__2; ll++) aa[ll] *= t;
40     kp1 = k + 1;
41     if (n < kp1) continue;
42     ax = aa;
43     for (j = kp1; j <= n; ++j) {
44       jn1        = j*n;
45       t          = a[k + jn1];
46       a[k + jn1] = 0.;
47       ay         = &a[1 + jn1];
48       for (ll=0; ll<k; ll++) ay[ll] += t*ax[ll];
49     }
50   }
51 
52   /*    form inverse(u)*inverse(l) */
53 
54   nm1 = n - 1;
55   if (nm1 < 1) {
56       PetscFunctionReturn(0);
57   }
58   for (kb = 1; kb <= nm1; ++kb) {
59     k   = n - kb;
60     kn  = k*n;
61     kp1 = k + 1;
62     aa  = a + kn;
63     for (i = kp1; i <= n; ++i) {
64       work[i] = aa[i];
65       aa[i]   = 0.;
66     }
67     for (j = kp1; j <= n; ++j) {
68       t  = work[j];
69       ax = &a[j * n + 1];
70       ay = &a[kn + 1];
71       for (ll=0; ll<n; ll++) ay[ll] += t*ax[ll];
72     }
73     l = ipvt[k];
74     if (l != k) {
75       ax = &a[kn + 1];
76       ay = &a[l * n + 1];
77       for (ll=0; ll<n; ll++) {
78         tmp    = ax[ll];
79         ax[ll] = ay[ll];
80         ay[ll] = tmp;
81       }
82     }
83   }
84   PetscFunctionReturn(0);
85 }
86 
87