Actual source code: baijsolvtran1.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <../src/mat/impls/baij/seq/baij.h>
  2:  #include <petsc/private/kernels/blockinvert.h>

  4: PetscErrorCode MatSolveTranspose_SeqBAIJ_1(Mat A,Vec bb,Vec xx)
  5: {
  6:   Mat_SeqBAIJ       *a    = (Mat_SeqBAIJ*)A->data;
  7:   IS                iscol = a->col,isrow = a->row;
  8:   PetscErrorCode    ierr;
  9:   const PetscInt    *rout,*cout,*r,*c,*adiag = a->diag,*ai = a->i,*aj = a->j,*vi;
 10:   PetscInt          i,n = a->mbs,j;
 11:   PetscInt          nz;
 12:   PetscScalar       *x,*tmp,s1;
 13:   const MatScalar   *aa = a->a,*v;
 14:   const PetscScalar *b;

 17:   VecGetArrayRead(bb,&b);
 18:   VecGetArray(xx,&x);
 19:   tmp  = a->solve_work;

 21:   ISGetIndices(isrow,&rout); r = rout;
 22:   ISGetIndices(iscol,&cout); c = cout;

 24:   /* copy the b into temp work space according to permutation */
 25:   for (i=0; i<n; i++) tmp[i] = b[c[i]];

 27:   /* forward solve the U^T */
 28:   for (i=0; i<n; i++) {
 29:     v   = aa + adiag[i+1] + 1;
 30:     vi  = aj + adiag[i+1] + 1;
 31:     nz  = adiag[i] - adiag[i+1] - 1;
 32:     s1  = tmp[i];
 33:     s1 *= v[nz];  /* multiply by inverse of diagonal entry */
 34:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
 35:     tmp[i] = s1;
 36:   }

 38:   /* backward solve the L^T */
 39:   for (i=n-1; i>=0; i--) {
 40:     v  = aa + ai[i];
 41:     vi = aj + ai[i];
 42:     nz = ai[i+1] - ai[i];
 43:     s1 = tmp[i];
 44:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
 45:   }

 47:   /* copy tmp into x according to permutation */
 48:   for (i=0; i<n; i++) x[r[i]] = tmp[i];

 50:   ISRestoreIndices(isrow,&rout);
 51:   ISRestoreIndices(iscol,&cout);
 52:   VecRestoreArrayRead(bb,&b);
 53:   VecRestoreArray(xx,&x);

 55:   PetscLogFlops(2.0*a->nz-A->cmap->n);
 56:   return(0);
 57: }

 59: PetscErrorCode MatSolveTranspose_SeqBAIJ_1_inplace(Mat A,Vec bb,Vec xx)
 60: {
 61:   Mat_SeqBAIJ       *a   =(Mat_SeqBAIJ*)A->data;
 62:   IS                iscol=a->col,isrow=a->row;
 63:   PetscErrorCode    ierr;
 64:   const PetscInt    *r,*c,*rout,*cout;
 65:   const PetscInt    *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
 66:   PetscInt          i,nz;
 67:   const MatScalar   *aa=a->a,*v;
 68:   PetscScalar       s1,*x,*t;
 69:   const PetscScalar *b;

 72:   VecGetArrayRead(bb,&b);
 73:   VecGetArray(xx,&x);
 74:   t    = a->solve_work;

 76:   ISGetIndices(isrow,&rout); r = rout;
 77:   ISGetIndices(iscol,&cout); c = cout;

 79:   /* copy the b into temp work space according to permutation */
 80:   for (i=0; i<n; i++) t[i] = b[c[i]];

 82:   /* forward solve the U^T */
 83:   for (i=0; i<n; i++) {

 85:     v = aa + diag[i];
 86:     /* multiply by the inverse of the block diagonal */
 87:     s1 = (*v++)*t[i];
 88:     vi = aj + diag[i] + 1;
 89:     nz = ai[i+1] - diag[i] - 1;
 90:     while (nz--) {
 91:       t[*vi++] -= (*v++)*s1;
 92:     }
 93:     t[i] = s1;
 94:   }
 95:   /* backward solve the L^T */
 96:   for (i=n-1; i>=0; i--) {
 97:     v  = aa + diag[i] - 1;
 98:     vi = aj + diag[i] - 1;
 99:     nz = diag[i] - ai[i];
100:     s1 = t[i];
101:     while (nz--) {
102:       t[*vi--] -=  (*v--)*s1;
103:     }
104:   }

106:   /* copy t into x according to permutation */
107:   for (i=0; i<n; i++) x[r[i]] = t[i];

109:   ISRestoreIndices(isrow,&rout);
110:   ISRestoreIndices(iscol,&cout);
111:   VecRestoreArrayRead(bb,&b);
112:   VecRestoreArray(xx,&x);
113:   PetscLogFlops(2.0*(a->nz) - A->cmap->n);
114:   return(0);
115: }