xref: /petsc/src/mat/impls/baij/seq/baijfact4.c (revision e8c0849ab8fe171bed529bea27238c9b402db591)
183287d42SBarry Smith /*
283287d42SBarry Smith     Factorization code for BAIJ format.
383287d42SBarry Smith */
4c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
5af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
683287d42SBarry Smith 
MatILUFactorNumeric_SeqBAIJ_N_inplace(Mat C,Mat A,const MatFactorInfo * info)7*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_N_inplace(Mat C, Mat A, const MatFactorInfo *info)
8d71ae5a4SJacob Faibussowitsch {
983287d42SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
1083287d42SBarry Smith   IS              isrow = b->row, isicol = b->icol;
115d0c19d7SBarry Smith   const PetscInt *r, *ic;
125d0c19d7SBarry Smith   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
133bc0b13bSBarry Smith   PetscInt       *ajtmpold, *ajtmp, nz, row, *ai = a->i, *aj = a->j, k, flg;
14*421480d9SBarry Smith   const PetscInt *diag_offset;
15*421480d9SBarry Smith   PetscInt        diag, bs = A->rmap->bs, bs2 = a->bs2, *pj, *v_pivots;
1683287d42SBarry Smith   MatScalar      *ba = b->a, *aa = a->a, *pv, *v, *rtmp, *multiplier, *v_work, *pc, *w;
175f8bbccaSHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
1883287d42SBarry Smith 
1983287d42SBarry Smith   PetscFunctionBegin;
20*421480d9SBarry Smith   /* Since A is C and C is labeled as a factored matrix we need to lie to MatGetDiagonalMarkers_SeqBAIJ() to get it to compute the diagonals */
21*421480d9SBarry Smith   A->factortype = MAT_FACTOR_NONE;
22*421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &diag_offset, NULL));
23*421480d9SBarry Smith   A->factortype = MAT_FACTOR_ILU;
249566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
259566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
265f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
275f8bbccaSHong Zhang 
289566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(bs2 * (n + 1), &rtmp));
2983287d42SBarry Smith   /* generate work space needed by dense LU factorization */
309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs, &v_work, bs2, &multiplier, bs, &v_pivots));
3183287d42SBarry Smith 
3283287d42SBarry Smith   for (i = 0; i < n; i++) {
3383287d42SBarry Smith     nz    = bi[i + 1] - bi[i];
3483287d42SBarry Smith     ajtmp = bj + bi[i];
3548a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * ajtmp[j], bs2));
3683287d42SBarry Smith     /* load in initial (unfactored row) */
3783287d42SBarry Smith     nz       = ai[r[i] + 1] - ai[r[i]];
3883287d42SBarry Smith     ajtmpold = aj + ai[r[i]];
3983287d42SBarry Smith     v        = aa + bs2 * ai[r[i]];
4048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmpold[j]], v + bs2 * j, bs2));
4183287d42SBarry Smith     row = *ajtmp++;
4283287d42SBarry Smith     while (row < i) {
4383287d42SBarry Smith       pc = rtmp + bs2 * row;
4483287d42SBarry Smith       /*      if (*pc) { */
45c35f09e5SBarry Smith       for (flg = 0, k = 0; k < bs2; k++) {
46c35f09e5SBarry Smith         if (pc[k] != 0.0) {
47c35f09e5SBarry Smith           flg = 1;
48c35f09e5SBarry Smith           break;
49c35f09e5SBarry Smith         }
50c35f09e5SBarry Smith       }
5183287d42SBarry Smith       if (flg) {
5283287d42SBarry Smith         pv = ba + bs2 * diag_offset[row];
5383287d42SBarry Smith         pj = bj + diag_offset[row] + 1;
5496b95a6bSBarry Smith         PetscKernel_A_gets_A_times_B(bs, pc, pv, multiplier);
5583287d42SBarry Smith         nz = bi[row + 1] - diag_offset[row] - 1;
5683287d42SBarry Smith         pv += bs2;
57ad540459SPierre Jolivet         for (j = 0; j < nz; j++) PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j);
589566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * bs * bs2 * (nz + 1.0) - bs));
5983287d42SBarry Smith       }
6083287d42SBarry Smith       row = *ajtmp++;
6183287d42SBarry Smith     }
6283287d42SBarry Smith     /* finished row so stick it into b->a */
6383287d42SBarry Smith     pv = ba + bs2 * bi[i];
6483287d42SBarry Smith     pj = bj + bi[i];
6583287d42SBarry Smith     nz = bi[i + 1] - bi[i];
6648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
6783287d42SBarry Smith     diag = diag_offset[i] - bi[i];
6883287d42SBarry Smith     /* invert diagonal block */
6983287d42SBarry Smith     w = pv + bs2 * diag;
705f8bbccaSHong Zhang 
719566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs, w, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
727b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
7383287d42SBarry Smith   }
7483287d42SBarry Smith 
759566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
769566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v_work, multiplier, v_pivots));
779566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
789566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
7926fbe8dcSKarl Rupp 
8006e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_N_inplace;
8106e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N_inplace;
8283287d42SBarry Smith   C->assembled           = PETSC_TRUE;
8326fbe8dcSKarl Rupp 
849566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8683287d42SBarry Smith }
87