xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
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>
172bb67c89SBarry Smith 
184a2ae208SSatish Balay #undef __FUNCT__
1996b95a6bSBarry Smith #define __FUNCT__ "PetscLINPACKgedi"
2096b95a6bSBarry Smith PetscErrorCode PetscLINPACKgedi(MatScalar *a,PetscInt n,PetscInt *ipvt,MatScalar *work)
212bb67c89SBarry Smith {
22690b6cddSBarry Smith   PetscInt  i__2,kb,kp1,nm1,i,j,k,l,ll,kn,knp1,jn1;
233f1db9ecSBarry Smith   MatScalar *aa,*ax,*ay,tmp;
243f1db9ecSBarry Smith   MatScalar t;
252bb67c89SBarry Smith 
263a40ed3dSBarry Smith   PetscFunctionBegin;
276d84be18SBarry Smith   --work;
286d84be18SBarry Smith   --ipvt;
2939d66777SBarry Smith   a -= n + 1;
302bb67c89SBarry Smith 
312bb67c89SBarry Smith   /*     compute inverse(u) */
322bb67c89SBarry Smith 
336d84be18SBarry Smith   for (k = 1; k <= n; ++k) {
3439d66777SBarry Smith     kn      = k*n;
3539d66777SBarry Smith     knp1    = kn + k;
3639d66777SBarry Smith     a[knp1] = 1.0 / a[knp1];
3739d66777SBarry Smith     t       = -a[knp1];
382bb67c89SBarry Smith     i__2    = k - 1;
3939d66777SBarry Smith     aa      = &a[1 + kn];
406d84be18SBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= t;
412bb67c89SBarry Smith     kp1 = k + 1;
426d84be18SBarry Smith     if (n < kp1) continue;
436d84be18SBarry Smith     ax = aa;
446d84be18SBarry Smith     for (j = kp1; j <= n; ++j) {
458d3e6ddaSBarry Smith       jn1        = j*n;
468d3e6ddaSBarry Smith       t          = a[k + jn1];
478d3e6ddaSBarry Smith       a[k + jn1] = 0.;
488d3e6ddaSBarry Smith       ay         = &a[1 + jn1];
49*26fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += t*ax[ll];
502bb67c89SBarry Smith     }
512bb67c89SBarry Smith   }
522bb67c89SBarry Smith 
532bb67c89SBarry Smith   /*    form inverse(u)*inverse(l) */
542bb67c89SBarry Smith 
556d84be18SBarry Smith   nm1 = n - 1;
562bb67c89SBarry Smith   if (nm1 < 1) {
573a40ed3dSBarry Smith       PetscFunctionReturn(0);
582bb67c89SBarry Smith   }
596d84be18SBarry Smith   for (kb = 1; kb <= nm1; ++kb) {
606d84be18SBarry Smith     k   = n - kb;
6139d66777SBarry Smith     kn  = k*n;
622bb67c89SBarry Smith     kp1 = k + 1;
6339d66777SBarry Smith     aa  = a + kn;
646d84be18SBarry Smith     for (i = kp1; i <= n; ++i) {
656d84be18SBarry Smith       work[i] = aa[i];
666d84be18SBarry Smith       aa[i]   = 0.;
672bb67c89SBarry Smith     }
686d84be18SBarry Smith     for (j = kp1; j <= n; ++j) {
692bb67c89SBarry Smith       t  = work[j];
706d84be18SBarry Smith       ax = &a[j * n + 1];
7139d66777SBarry Smith       ay = &a[kn + 1];
72*26fbe8dcSKarl Rupp       for (ll=0; ll<n; ll++) ay[ll] += t*ax[ll];
732bb67c89SBarry Smith     }
742bb67c89SBarry Smith     l = ipvt[k];
752bb67c89SBarry Smith     if (l != k) {
7639d66777SBarry Smith       ax = &a[kn + 1];
776d84be18SBarry Smith       ay = &a[l * n + 1];
786d84be18SBarry Smith       for (ll=0; ll<n; ll++) {
796d84be18SBarry Smith         tmp    = ax[ll];
806d84be18SBarry Smith         ax[ll] = ay[ll];
816d84be18SBarry Smith         ay[ll] = tmp;
822bb67c89SBarry Smith       }
832bb67c89SBarry Smith     }
846d84be18SBarry Smith   }
853a40ed3dSBarry Smith   PetscFunctionReturn(0);
866d84be18SBarry Smith }
872bb67c89SBarry Smith 
88