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