Actual source code: baijfact4.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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: /* ----------------------------------------------------------- */
 11: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N_inplace(Mat C,Mat A,const MatFactorInfo *info)
 12: {
 13:   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
 14:   IS             isrow = b->row,isicol = b->icol;
 16:   const PetscInt *r,*ic;
 17:   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
 18:   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j,k,flg;
 19:   PetscInt       *diag_offset=b->diag,diag,bs=A->rmap->bs,bs2 = a->bs2,*pj,*v_pivots;
 20:   MatScalar      *ba = b->a,*aa = a->a,*pv,*v,*rtmp,*multiplier,*v_work,*pc,*w;
 21:   PetscBool      allowzeropivot,zeropivotdetected;

 24:   ISGetIndices(isrow,&r);
 25:   ISGetIndices(isicol,&ic);
 26:   allowzeropivot = PetscNot(A->erroriffailure);

 28:   PetscMalloc1(bs2*(n+1),&rtmp);
 29:   PetscMemzero(rtmp,(bs2*n+1)*sizeof(MatScalar));
 30:   /* generate work space needed by dense LU factorization */
 31:   PetscMalloc3(bs,&v_work,bs2,&multiplier,bs,&v_pivots);

 33:   for (i=0; i<n; i++) {
 34:     nz    = bi[i+1] - bi[i];
 35:     ajtmp = bj + bi[i];
 36:     for  (j=0; j<nz; j++) {
 37:       PetscMemzero(rtmp+bs2*ajtmp[j],bs2*sizeof(MatScalar));
 38:     }
 39:     /* load in initial (unfactored row) */
 40:     nz       = ai[r[i]+1] - ai[r[i]];
 41:     ajtmpold = aj + ai[r[i]];
 42:     v        = aa + bs2*ai[r[i]];
 43:     for (j=0; j<nz; j++) {
 44:       PetscMemcpy(rtmp+bs2*ic[ajtmpold[j]],v+bs2*j,bs2*sizeof(MatScalar));
 45:     }
 46:     row = *ajtmp++;
 47:     while (row < i) {
 48:       pc = rtmp + bs2*row;
 49: /*      if (*pc) { */
 50:       for (flg=0,k=0; k<bs2; k++) {
 51:         if (pc[k]!=0.0) {
 52:           flg = 1;
 53:           break;
 54:         }
 55:       }
 56:       if (flg) {
 57:         pv = ba + bs2*diag_offset[row];
 58:         pj = bj + diag_offset[row] + 1;
 59:         PetscKernel_A_gets_A_times_B(bs,pc,pv,multiplier);
 60:         nz  = bi[row+1] - diag_offset[row] - 1;
 61:         pv += bs2;
 62:         for (j=0; j<nz; j++) {
 63:           PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
 64:         }
 65:         PetscLogFlops(2.0*bs*bs2*(nz+1.0)-bs);
 66:       }
 67:       row = *ajtmp++;
 68:     }
 69:     /* finished row so stick it into b->a */
 70:     pv = ba + bs2*bi[i];
 71:     pj = bj + bi[i];
 72:     nz = bi[i+1] - bi[i];
 73:     for (j=0; j<nz; j++) {
 74:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
 75:     }
 76:     diag = diag_offset[i] - bi[i];
 77:     /* invert diagonal block */
 78:     w    = pv + bs2*diag;

 80:     PetscKernel_A_gets_inverse_A(bs,w,v_pivots,v_work,allowzeropivot,&zeropivotdetected);
 81:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
 82:   }

 84:   PetscFree(rtmp);
 85:   PetscFree3(v_work,multiplier,v_pivots);
 86:   ISRestoreIndices(isicol,&ic);
 87:   ISRestoreIndices(isrow,&r);

 89:   C->ops->solve          = MatSolve_SeqBAIJ_N_inplace;
 90:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N_inplace;
 91:   C->assembled           = PETSC_TRUE;

 93:   PetscLogFlops(1.333333333333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
 94:   return(0);
 95: }