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