xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision ad540459ab38c4a232710a68077e487eb6627fe2)
1be1d678aSKris Buschelman 
26d84be18SBarry Smith /*
36d84be18SBarry Smith               This file creating by running f2c
46d84be18SBarry Smith             linpack. this version dated 08/14/78
56d84be18SBarry Smith       cleve moler, university of new mexico, argonne national lab.
66d84be18SBarry Smith 
76d84be18SBarry Smith       Computes the inverse of a matrix given its factors and pivots
896b95a6bSBarry Smith     calculated by PetscLINPACKgefa(). Performed in-place for an n by n
99d8b60e5SBarry Smith     dense matrix.
109d8b60e5SBarry Smith 
119d8b60e5SBarry Smith        Used by the sparse factorization routines in
12dd882469SBarry Smith      src/mat/impls/baij/seq
1371c5468dSBarry Smith 
142bb67c89SBarry Smith */
152bb67c89SBarry Smith 
16c6db04a5SJed Brown #include <petscsys.h>
175f748d31SJed Brown #include <petsc/private/kernels/blockinvert.h>
182bb67c89SBarry Smith 
199371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscLINPACKgedi(MatScalar *a, PetscInt n, PetscInt *ipvt, MatScalar *work) {
20690b6cddSBarry Smith   PetscInt   i__2, kb, kp1, nm1, i, j, k, l, ll, kn, knp1, jn1;
213f1db9ecSBarry Smith   MatScalar *aa, *ax, *ay, tmp;
223f1db9ecSBarry Smith   MatScalar  t;
232bb67c89SBarry Smith 
243a40ed3dSBarry Smith   PetscFunctionBegin;
256d84be18SBarry Smith   --work;
266d84be18SBarry Smith   --ipvt;
2739d66777SBarry Smith   a -= n + 1;
282bb67c89SBarry Smith 
292bb67c89SBarry Smith   /*     compute inverse(u) */
302bb67c89SBarry Smith 
316d84be18SBarry Smith   for (k = 1; k <= n; ++k) {
3239d66777SBarry Smith     kn      = k * n;
3339d66777SBarry Smith     knp1    = kn + k;
3439d66777SBarry Smith     a[knp1] = 1.0 / a[knp1];
3539d66777SBarry Smith     t       = -a[knp1];
362bb67c89SBarry Smith     i__2    = k - 1;
3739d66777SBarry Smith     aa      = &a[1 + kn];
386d84be18SBarry Smith     for (ll = 0; ll < i__2; ll++) aa[ll] *= t;
392bb67c89SBarry Smith     kp1 = k + 1;
406d84be18SBarry Smith     if (n < kp1) continue;
416d84be18SBarry Smith     ax = aa;
426d84be18SBarry Smith     for (j = kp1; j <= n; ++j) {
438d3e6ddaSBarry Smith       jn1        = j * n;
448d3e6ddaSBarry Smith       t          = a[k + jn1];
458d3e6ddaSBarry Smith       a[k + jn1] = 0.;
468d3e6ddaSBarry Smith       ay         = &a[1 + jn1];
4726fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += t * ax[ll];
482bb67c89SBarry Smith     }
492bb67c89SBarry Smith   }
502bb67c89SBarry Smith 
512bb67c89SBarry Smith   /*    form inverse(u)*inverse(l) */
522bb67c89SBarry Smith 
536d84be18SBarry Smith   nm1 = n - 1;
54*ad540459SPierre Jolivet   if (nm1 < 1) PetscFunctionReturn(0);
556d84be18SBarry Smith   for (kb = 1; kb <= nm1; ++kb) {
566d84be18SBarry Smith     k   = n - kb;
5739d66777SBarry Smith     kn  = k * n;
582bb67c89SBarry Smith     kp1 = k + 1;
5939d66777SBarry Smith     aa  = a + kn;
606d84be18SBarry Smith     for (i = kp1; i <= n; ++i) {
616d84be18SBarry Smith       work[i] = aa[i];
626d84be18SBarry Smith       aa[i]   = 0.;
632bb67c89SBarry Smith     }
646d84be18SBarry Smith     for (j = kp1; j <= n; ++j) {
652bb67c89SBarry Smith       t  = work[j];
666d84be18SBarry Smith       ax = &a[j * n + 1];
6739d66777SBarry Smith       ay = &a[kn + 1];
6826fbe8dcSKarl Rupp       for (ll = 0; ll < n; ll++) ay[ll] += t * ax[ll];
692bb67c89SBarry Smith     }
702bb67c89SBarry Smith     l = ipvt[k];
712bb67c89SBarry Smith     if (l != k) {
7239d66777SBarry Smith       ax = &a[kn + 1];
736d84be18SBarry Smith       ay = &a[l * n + 1];
746d84be18SBarry Smith       for (ll = 0; ll < n; ll++) {
756d84be18SBarry Smith         tmp    = ax[ll];
766d84be18SBarry Smith         ax[ll] = ay[ll];
776d84be18SBarry Smith         ay[ll] = tmp;
782bb67c89SBarry Smith       }
792bb67c89SBarry Smith     }
806d84be18SBarry Smith   }
813a40ed3dSBarry Smith   PetscFunctionReturn(0);
826d84be18SBarry Smith }
83