Actual source code: dgedi.c
petsc-3.7.3 2016-08-01
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.
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.
11: Used by the sparse factorization routines in
12: src/mat/impls/baij/seq
14: */
16: #include <petscsys.h>
20: PETSC_INTERN 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;
27: --work;
28: --ipvt;
29: a -= n + 1;
31: /* compute inverse(u) */
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++) ay[ll] += t*ax[ll];
50: }
51: }
53: /* form inverse(u)*inverse(l) */
55: nm1 = n - 1;
56: if (nm1 < 1) {
57: return(0);
58: }
59: for (kb = 1; kb <= nm1; ++kb) {
60: k = n - kb;
61: kn = k*n;
62: kp1 = k + 1;
63: aa = a + kn;
64: for (i = kp1; i <= n; ++i) {
65: work[i] = aa[i];
66: aa[i] = 0.;
67: }
68: for (j = kp1; j <= n; ++j) {
69: t = work[j];
70: ax = &a[j * n + 1];
71: ay = &a[kn + 1];
72: for (ll=0; ll<n; ll++) ay[ll] += t*ax[ll];
73: }
74: l = ipvt[k];
75: if (l != k) {
76: ax = &a[kn + 1];
77: ay = &a[l * n + 1];
78: for (ll=0; ll<n; ll++) {
79: tmp = ax[ll];
80: ax[ll] = ay[ll];
81: ay[ll] = tmp;
82: }
83: }
84: }
85: return(0);
86: }