1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*71c5468dSBarry Smith static char vcid[] = "$Id: dgedi.c,v 1.10 1998/12/24 02:56:37 bsmith Exp bsmith $"; 36d84be18SBarry Smith #endif 46d84be18SBarry Smith 56d84be18SBarry Smith /* 66d84be18SBarry Smith This file creating by running f2c 76d84be18SBarry Smith linpack. this version dated 08/14/78 86d84be18SBarry Smith cleve moler, university of new mexico, argonne national lab. 96d84be18SBarry Smith 106d84be18SBarry Smith Computes the inverse of a matrix given its factors and pivots 119d8b60e5SBarry Smith calculated by Linpack_DGEFA(). Performed in-place for an n by n 129d8b60e5SBarry Smith dense matrix. 139d8b60e5SBarry Smith 149d8b60e5SBarry Smith Used by the sparse factorization routines in 159d8b60e5SBarry Smith src/mat/impls/baij/seq and src/mat/impls/bdiag/seq 16*71c5468dSBarry Smith 17*71c5468dSBarry Smith See also src/inline/ilu.h 182bb67c89SBarry Smith */ 192bb67c89SBarry Smith 206d84be18SBarry Smith #include "petsc.h" 212bb67c89SBarry Smith 225615d1e5SSatish Balay #undef __FUNC__ 235615d1e5SSatish Balay #define __FUNC__ "Linpack_DGEDI" 243f1db9ecSBarry Smith int Linpack_DGEDI(MatScalar *a,int n,int *ipvt,MatScalar *work) 252bb67c89SBarry Smith { 268d3e6ddaSBarry Smith int i__2,kb, kp1, nm1,i, j, k, l, ll,kn,knp1,jn1; 273f1db9ecSBarry Smith MatScalar *aa,*ax,*ay,tmp; 283f1db9ecSBarry Smith MatScalar t; 292bb67c89SBarry Smith 303a40ed3dSBarry Smith PetscFunctionBegin; 316d84be18SBarry Smith --work; 326d84be18SBarry Smith --ipvt; 3339d66777SBarry Smith a -= n + 1; 342bb67c89SBarry Smith 352bb67c89SBarry Smith /* compute inverse(u) */ 362bb67c89SBarry Smith 376d84be18SBarry Smith for (k = 1; k <= n; ++k) { 3839d66777SBarry Smith kn = k*n; 3939d66777SBarry Smith knp1 = kn + k; 4039d66777SBarry Smith a[knp1] = 1.0 / a[knp1]; 4139d66777SBarry Smith t = -a[knp1]; 422bb67c89SBarry Smith i__2 = k - 1; 4339d66777SBarry Smith aa = &a[1 + kn]; 446d84be18SBarry Smith for ( ll=0; ll<i__2; ll++ ) aa[ll] *= t; 452bb67c89SBarry Smith kp1 = k + 1; 466d84be18SBarry Smith if (n < kp1) continue; 476d84be18SBarry Smith ax = aa; 486d84be18SBarry Smith for (j = kp1; j <= n; ++j) { 498d3e6ddaSBarry Smith jn1 = j*n; 508d3e6ddaSBarry Smith t = a[k + jn1]; 518d3e6ddaSBarry Smith a[k + jn1] = 0.; 528d3e6ddaSBarry Smith ay = &a[1 + jn1]; 536d84be18SBarry Smith for ( ll=0; ll<k; ll++ ) { 546d84be18SBarry Smith ay[ll] += t*ax[ll]; 552bb67c89SBarry Smith } 562bb67c89SBarry Smith } 572bb67c89SBarry Smith } 582bb67c89SBarry Smith 592bb67c89SBarry Smith /* form inverse(u)*inverse(l) */ 602bb67c89SBarry Smith 616d84be18SBarry Smith nm1 = n - 1; 622bb67c89SBarry Smith if (nm1 < 1) { 633a40ed3dSBarry Smith PetscFunctionReturn(0); 642bb67c89SBarry Smith } 656d84be18SBarry Smith for (kb = 1; kb <= nm1; ++kb) { 666d84be18SBarry Smith k = n - kb; 6739d66777SBarry Smith kn = k*n; 682bb67c89SBarry Smith kp1 = k + 1; 6939d66777SBarry Smith aa = a + kn; 706d84be18SBarry Smith for (i = kp1; i <= n; ++i) { 716d84be18SBarry Smith work[i] = aa[i]; 726d84be18SBarry Smith aa[i] = 0.; 732bb67c89SBarry Smith } 746d84be18SBarry Smith for (j = kp1; j <= n; ++j) { 752bb67c89SBarry Smith t = work[j]; 766d84be18SBarry Smith ax = &a[j * n + 1]; 7739d66777SBarry Smith ay = &a[kn + 1]; 786d84be18SBarry Smith for ( ll=0; ll<n; ll++ ) { 796d84be18SBarry Smith ay[ll] += t*ax[ll]; 806d84be18SBarry Smith } 812bb67c89SBarry Smith } 822bb67c89SBarry Smith l = ipvt[k]; 832bb67c89SBarry Smith if (l != k) { 8439d66777SBarry Smith ax = &a[kn + 1]; 856d84be18SBarry Smith ay = &a[l * n + 1]; 866d84be18SBarry Smith for ( ll=0; ll<n; ll++ ) { 876d84be18SBarry Smith tmp = ax[ll]; 886d84be18SBarry Smith ax[ll] = ay[ll]; 896d84be18SBarry Smith ay[ll] = tmp; 902bb67c89SBarry Smith } 912bb67c89SBarry Smith } 926d84be18SBarry Smith } 933a40ed3dSBarry Smith PetscFunctionReturn(0); 946d84be18SBarry Smith } 952bb67c89SBarry Smith 96