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