xref: /petsc/src/mat/impls/baij/seq/baijfact4.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 /*
2     Factorization code for BAIJ format.
3 */
4 #include <../src/mat/impls/baij/seq/baij.h>
5 #include <petsc/private/kernels/blockinvert.h>
6 
7 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N_inplace(Mat C, Mat A, const MatFactorInfo *info)
8 {
9   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
10   IS              isrow = b->row, isicol = b->icol;
11   const PetscInt *r, *ic;
12   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
13   PetscInt       *ajtmpold, *ajtmp, nz, row, *ai = a->i, *aj = a->j, k, flg;
14   PetscInt       *diag_offset = b->diag, diag, bs = A->rmap->bs, bs2 = a->bs2, *pj, *v_pivots;
15   MatScalar      *ba = b->a, *aa = a->a, *pv, *v, *rtmp, *multiplier, *v_work, *pc, *w;
16   PetscBool       allowzeropivot, zeropivotdetected;
17 
18   PetscFunctionBegin;
19   PetscCall(ISGetIndices(isrow, &r));
20   PetscCall(ISGetIndices(isicol, &ic));
21   allowzeropivot = PetscNot(A->erroriffailure);
22 
23   PetscCall(PetscCalloc1(bs2 * (n + 1), &rtmp));
24   /* generate work space needed by dense LU factorization */
25   PetscCall(PetscMalloc3(bs, &v_work, bs2, &multiplier, bs, &v_pivots));
26 
27   for (i = 0; i < n; i++) {
28     nz    = bi[i + 1] - bi[i];
29     ajtmp = bj + bi[i];
30     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * ajtmp[j], bs2));
31     /* load in initial (unfactored row) */
32     nz       = ai[r[i] + 1] - ai[r[i]];
33     ajtmpold = aj + ai[r[i]];
34     v        = aa + bs2 * ai[r[i]];
35     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmpold[j]], v + bs2 * j, bs2));
36     row = *ajtmp++;
37     while (row < i) {
38       pc = rtmp + bs2 * row;
39       /*      if (*pc) { */
40       for (flg = 0, k = 0; k < bs2; k++) {
41         if (pc[k] != 0.0) {
42           flg = 1;
43           break;
44         }
45       }
46       if (flg) {
47         pv = ba + bs2 * diag_offset[row];
48         pj = bj + diag_offset[row] + 1;
49         PetscKernel_A_gets_A_times_B(bs, pc, pv, multiplier);
50         nz = bi[row + 1] - diag_offset[row] - 1;
51         pv += bs2;
52         for (j = 0; j < nz; j++) PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j);
53         PetscCall(PetscLogFlops(2.0 * bs * bs2 * (nz + 1.0) - bs));
54       }
55       row = *ajtmp++;
56     }
57     /* finished row so stick it into b->a */
58     pv = ba + bs2 * bi[i];
59     pj = bj + bi[i];
60     nz = bi[i + 1] - bi[i];
61     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
62     diag = diag_offset[i] - bi[i];
63     /* invert diagonal block */
64     w = pv + bs2 * diag;
65 
66     PetscCall(PetscKernel_A_gets_inverse_A(bs, w, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
67     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
68   }
69 
70   PetscCall(PetscFree(rtmp));
71   PetscCall(PetscFree3(v_work, multiplier, v_pivots));
72   PetscCall(ISRestoreIndices(isicol, &ic));
73   PetscCall(ISRestoreIndices(isrow, &r));
74 
75   C->ops->solve          = MatSolve_SeqBAIJ_N_inplace;
76   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N_inplace;
77   C->assembled           = PETSC_TRUE;
78 
79   PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82