Actual source code: baijsolvtrannat1.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <../src/mat/impls/baij/seq/baij.h>

  3: PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)
  4: {
  5:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
  6:   PetscErrorCode    ierr;
  7:   const PetscInt    *adiag = a->diag,*ai = a->i,*aj = a->j,*vi;
  8:   PetscInt          i,n = a->mbs,j;
  9:   PetscInt          nz;
 10:   PetscScalar       *x,*tmp,s1;
 11:   const MatScalar   *aa = a->a,*v;
 12:   const PetscScalar *b;

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


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

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

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

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

 46:   VecRestoreArrayRead(bb,&b);
 47:   VecRestoreArray(xx,&x);

 49:   PetscLogFlops(2.0*a->nz-A->cmap->n);
 50:   return(0);
 51: }

 53: PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
 54: {
 55:   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data;
 56:   PetscErrorCode  ierr;
 57:   PetscInt        i,nz;
 58:   const PetscInt  *diag = a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
 59:   const MatScalar *aa   =a->a,*v;
 60:   PetscScalar     s1,*x;

 63:   VecCopy(bb,xx);
 64:   VecGetArray(xx,&x);

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

 69:     v = aa + diag[i];
 70:     /* multiply by the inverse of the block diagonal */
 71:     s1 = (*v++)*x[i];
 72:     vi = aj + diag[i] + 1;
 73:     nz = ai[i+1] - diag[i] - 1;
 74:     while (nz--) {
 75:       x[*vi++] -= (*v++)*s1;
 76:     }
 77:     x[i] = s1;
 78:   }
 79:   /* backward solve the L^T */
 80:   for (i=n-1; i>=0; i--) {
 81:     v  = aa + diag[i] - 1;
 82:     vi = aj + diag[i] - 1;
 83:     nz = diag[i] - ai[i];
 84:     s1 = x[i];
 85:     while (nz--) {
 86:       x[*vi--] -=  (*v--)*s1;
 87:     }
 88:   }
 89:   VecRestoreArray(xx,&x);
 90:   PetscLogFlops(2.0*(a->nz) - A->cmap->n);
 91:   return(0);
 92: }