xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision 96b95a6bf9ed541503cb02ce9899069044d5ddfa)
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
8*96b95a6bSBarry 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__
19*96b95a6bSBarry Smith #define __FUNCT__ "PetscLINPACKgedi"
20*96b95a6bSBarry 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];
496d84be18SBarry Smith             for (ll=0; ll<k; ll++) {
506d84be18SBarry Smith               ay[ll] += t*ax[ll];
512bb67c89SBarry Smith             }
522bb67c89SBarry Smith 	}
532bb67c89SBarry Smith     }
542bb67c89SBarry Smith 
552bb67c89SBarry Smith    /*    form inverse(u)*inverse(l) */
562bb67c89SBarry Smith 
576d84be18SBarry Smith     nm1 = n - 1;
582bb67c89SBarry Smith     if (nm1 < 1) {
593a40ed3dSBarry Smith 	PetscFunctionReturn(0);
602bb67c89SBarry Smith     }
616d84be18SBarry Smith     for (kb = 1; kb <= nm1; ++kb) {
626d84be18SBarry Smith 	k   = n - kb;
6339d66777SBarry Smith         kn  = k*n;
642bb67c89SBarry Smith 	kp1 = k + 1;
6539d66777SBarry Smith         aa  = a + kn;
666d84be18SBarry Smith 	for (i = kp1; i <= n; ++i) {
676d84be18SBarry Smith 	    work[i] = aa[i];
686d84be18SBarry Smith 	    aa[i]   = 0.;
692bb67c89SBarry Smith 	}
706d84be18SBarry Smith 	for (j = kp1; j <= n; ++j) {
712bb67c89SBarry Smith 	    t = work[j];
726d84be18SBarry Smith             ax = &a[j * n + 1];
7339d66777SBarry Smith             ay = &a[kn + 1];
746d84be18SBarry Smith             for (ll=0; ll<n; ll++) {
756d84be18SBarry Smith               ay[ll] += t*ax[ll];
766d84be18SBarry Smith             }
772bb67c89SBarry Smith 	}
782bb67c89SBarry Smith 	l = ipvt[k];
792bb67c89SBarry Smith 	if (l != k) {
8039d66777SBarry Smith             ax = &a[kn + 1];
816d84be18SBarry Smith             ay = &a[l * n + 1];
826d84be18SBarry Smith             for (ll=0; ll<n; ll++) {
836d84be18SBarry Smith               tmp    = ax[ll];
846d84be18SBarry Smith               ax[ll] = ay[ll];
856d84be18SBarry Smith               ay[ll] = tmp;
862bb67c89SBarry Smith             }
872bb67c89SBarry Smith 	}
886d84be18SBarry Smith     }
893a40ed3dSBarry Smith     PetscFunctionReturn(0);
906d84be18SBarry Smith }
912bb67c89SBarry Smith 
92