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