xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision 9de2952e64edfa51c4ec5537730b8a834bdd3686)
1a97b7df5SSatish Balay /*
29fd41dd0SBarry Smith        This routine was converted by f2c from Linpack source
39fd41dd0SBarry Smith              linpack. this version dated 08/14/78
49fd41dd0SBarry Smith       cleve moler, university of new mexico, argonne national lab.
559539b86SBarry Smith 
659539b86SBarry Smith         Does an LU factorization with partial pivoting of a dense
79d8b60e5SBarry Smith      n by n matrix.
89d8b60e5SBarry Smith 
99d8b60e5SBarry Smith        Used by the sparse factorization routines in
10dd882469SBarry Smith      src/mat/impls/baij/seq
1159539b86SBarry Smith 
129fd41dd0SBarry Smith */
13*9de2952eSStefano Zampini #include <petsc/private/kernels/blockinvert.h>
149fd41dd0SBarry Smith 
15*9de2952eSStefano Zampini PetscErrorCode PetscLINPACKgefa(MatScalar *a, PetscInt n, PetscInt *ipvt, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
16d71ae5a4SJacob Faibussowitsch {
17690b6cddSBarry Smith   PetscInt  i__2, i__3, kp1, nm1, j, k, l, ll, kn, knp1, jn1;
183f1db9ecSBarry Smith   MatScalar t, *ax, *ay, *aa;
19329f5518SBarry Smith   MatReal   tmp, max;
209fd41dd0SBarry Smith 
213a40ed3dSBarry Smith   PetscFunctionBegin;
225f8bbccaSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
235f8bbccaSHong Zhang 
249fd41dd0SBarry Smith   /* Parameter adjustments */
259fd41dd0SBarry Smith   --ipvt;
2639d66777SBarry Smith   a -= n + 1;
279fd41dd0SBarry Smith 
289fd41dd0SBarry Smith   /* Function Body */
299fd41dd0SBarry Smith   nm1 = n - 1;
3039d66777SBarry Smith   for (k = 1; k <= nm1; ++k) {
319fd41dd0SBarry Smith     kp1  = k + 1;
3239d66777SBarry Smith     kn   = k * n;
3339d66777SBarry Smith     knp1 = k * n + k;
349fd41dd0SBarry Smith 
359fd41dd0SBarry Smith     /* find l = pivot index */
369fd41dd0SBarry Smith 
379fd41dd0SBarry Smith     i__2 = n - k + 1;
3839d66777SBarry Smith     aa   = &a[knp1];
399fd41dd0SBarry Smith     max  = PetscAbsScalar(aa[0]);
409fd41dd0SBarry Smith     l    = 1;
419fd41dd0SBarry Smith     for (ll = 1; ll < i__2; ll++) {
429fd41dd0SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
439371c9d4SSatish Balay       if (tmp > max) {
449371c9d4SSatish Balay         max = tmp;
459371c9d4SSatish Balay         l   = ll + 1;
469371c9d4SSatish Balay       }
479fd41dd0SBarry Smith     }
489fd41dd0SBarry Smith     l += k - 1;
499fd41dd0SBarry Smith     ipvt[k] = l;
509fd41dd0SBarry Smith 
515f8bbccaSHong Zhang     if (a[l + kn] == 0.0) {
525f8bbccaSHong Zhang       if (allowzeropivot) {
539566063dSJacob Faibussowitsch         PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
5470d8d27cSBarry Smith         if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
5598921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
565f8bbccaSHong Zhang     }
579fd41dd0SBarry Smith 
589fd41dd0SBarry Smith     /* interchange if necessary */
599fd41dd0SBarry Smith     if (l != k) {
6039d66777SBarry Smith       t         = a[l + kn];
6139d66777SBarry Smith       a[l + kn] = a[knp1];
6239d66777SBarry Smith       a[knp1]   = t;
639fd41dd0SBarry Smith     }
649fd41dd0SBarry Smith 
659fd41dd0SBarry Smith     /* compute multipliers */
6639d66777SBarry Smith     t    = -1. / a[knp1];
679fd41dd0SBarry Smith     i__2 = n - k;
6839d66777SBarry Smith     aa   = &a[1 + knp1];
6926fbe8dcSKarl Rupp     for (ll = 0; ll < i__2; ll++) aa[ll] *= t;
709fd41dd0SBarry Smith 
719fd41dd0SBarry Smith     /* row elimination with column indexing */
7239d66777SBarry Smith     ax = aa;
739fd41dd0SBarry Smith     for (j = kp1; j <= n; ++j) {
748d3e6ddaSBarry Smith       jn1 = j * n;
758d3e6ddaSBarry Smith       t   = a[l + jn1];
769fd41dd0SBarry Smith       if (l != k) {
778d3e6ddaSBarry Smith         a[l + jn1] = a[k + jn1];
788d3e6ddaSBarry Smith         a[k + jn1] = t;
799fd41dd0SBarry Smith       }
809fd41dd0SBarry Smith 
819fd41dd0SBarry Smith       i__3 = n - k;
828d3e6ddaSBarry Smith       ay   = &a[1 + k + jn1];
8326fbe8dcSKarl Rupp       for (ll = 0; ll < i__3; ll++) ay[ll] += t * ax[ll];
849fd41dd0SBarry Smith     }
859fd41dd0SBarry Smith   }
869fd41dd0SBarry Smith   ipvt[n] = n;
875f8bbccaSHong Zhang   if (a[n + n * n] == 0.0) {
885f8bbccaSHong Zhang     if (allowzeropivot) {
899566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", n - 1));
9070d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
9198921bdaSJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, n - 1);
925f8bbccaSHong Zhang   }
933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
949fd41dd0SBarry Smith }
95