Actual source code: dgedi.c
petsc-3.8.4 2018-03-24
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>
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;
25: --work;
26: --ipvt;
27: a -= n + 1;
29: /* compute inverse(u) */
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: }
51: /* form inverse(u)*inverse(l) */
53: nm1 = n - 1;
54: if (nm1 < 1) {
55: return(0);
56: }
57: for (kb = 1; kb <= nm1; ++kb) {
58: k = n - kb;
59: kn = k*n;
60: kp1 = k + 1;
61: aa = a + kn;
62: for (i = kp1; i <= n; ++i) {
63: work[i] = aa[i];
64: aa[i] = 0.;
65: }
66: for (j = kp1; j <= n; ++j) {
67: t = work[j];
68: ax = &a[j * n + 1];
69: ay = &a[kn + 1];
70: for (ll=0; ll<n; ll++) ay[ll] += t*ax[ll];
71: }
72: l = ipvt[k];
73: if (l != k) {
74: ax = &a[kn + 1];
75: ay = &a[l * n + 1];
76: for (ll=0; ll<n; ll++) {
77: tmp = ax[ll];
78: ax[ll] = ay[ll];
79: ay[ll] = tmp;
80: }
81: }
82: }
83: return(0);
84: }