1 #include <../src/mat/impls/baij/seq/baij.h> 2 #include <petsc/private/kernels/blockinvert.h> 3 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 const PetscInt *rout,*cout,*r,*c,*adiag = a->diag,*ai = a->i,*aj = a->j,*vi; 9 PetscInt i,n = a->mbs,j; 10 PetscInt nz; 11 PetscScalar *x,*tmp,s1; 12 const MatScalar *aa = a->a,*v; 13 const PetscScalar *b; 14 15 PetscFunctionBegin; 16 PetscCall(VecGetArrayRead(bb,&b)); 17 PetscCall(VecGetArray(xx,&x)); 18 tmp = a->solve_work; 19 20 PetscCall(ISGetIndices(isrow,&rout)); r = rout; 21 PetscCall(ISGetIndices(iscol,&cout)); c = cout; 22 23 /* copy the b into temp work space according to permutation */ 24 for (i=0; i<n; i++) tmp[i] = b[c[i]]; 25 26 /* forward solve the U^T */ 27 for (i=0; i<n; i++) { 28 v = aa + adiag[i+1] + 1; 29 vi = aj + adiag[i+1] + 1; 30 nz = adiag[i] - adiag[i+1] - 1; 31 s1 = tmp[i]; 32 s1 *= v[nz]; /* multiply by inverse of diagonal entry */ 33 for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j]; 34 tmp[i] = s1; 35 } 36 37 /* backward solve the L^T */ 38 for (i=n-1; i>=0; i--) { 39 v = aa + ai[i]; 40 vi = aj + ai[i]; 41 nz = ai[i+1] - ai[i]; 42 s1 = tmp[i]; 43 for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j]; 44 } 45 46 /* copy tmp into x according to permutation */ 47 for (i=0; i<n; i++) x[r[i]] = tmp[i]; 48 49 PetscCall(ISRestoreIndices(isrow,&rout)); 50 PetscCall(ISRestoreIndices(iscol,&cout)); 51 PetscCall(VecRestoreArrayRead(bb,&b)); 52 PetscCall(VecRestoreArray(xx,&x)); 53 54 PetscCall(PetscLogFlops(2.0*a->nz-A->cmap->n)); 55 PetscFunctionReturn(0); 56 } 57 58 PetscErrorCode MatSolveTranspose_SeqBAIJ_1_inplace(Mat A,Vec bb,Vec xx) 59 { 60 Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data; 61 IS iscol=a->col,isrow=a->row; 62 const PetscInt *r,*c,*rout,*cout; 63 const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j; 64 PetscInt i,nz; 65 const MatScalar *aa=a->a,*v; 66 PetscScalar s1,*x,*t; 67 const PetscScalar *b; 68 69 PetscFunctionBegin; 70 PetscCall(VecGetArrayRead(bb,&b)); 71 PetscCall(VecGetArray(xx,&x)); 72 t = a->solve_work; 73 74 PetscCall(ISGetIndices(isrow,&rout)); r = rout; 75 PetscCall(ISGetIndices(iscol,&cout)); c = cout; 76 77 /* copy the b into temp work space according to permutation */ 78 for (i=0; i<n; i++) t[i] = b[c[i]]; 79 80 /* forward solve the U^T */ 81 for (i=0; i<n; i++) { 82 83 v = aa + diag[i]; 84 /* multiply by the inverse of the block diagonal */ 85 s1 = (*v++)*t[i]; 86 vi = aj + diag[i] + 1; 87 nz = ai[i+1] - diag[i] - 1; 88 while (nz--) { 89 t[*vi++] -= (*v++)*s1; 90 } 91 t[i] = s1; 92 } 93 /* backward solve the L^T */ 94 for (i=n-1; i>=0; i--) { 95 v = aa + diag[i] - 1; 96 vi = aj + diag[i] - 1; 97 nz = diag[i] - ai[i]; 98 s1 = t[i]; 99 while (nz--) { 100 t[*vi--] -= (*v--)*s1; 101 } 102 } 103 104 /* copy t into x according to permutation */ 105 for (i=0; i<n; i++) x[r[i]] = t[i]; 106 107 PetscCall(ISRestoreIndices(isrow,&rout)); 108 PetscCall(ISRestoreIndices(iscol,&cout)); 109 PetscCall(VecRestoreArrayRead(bb,&b)); 110 PetscCall(VecRestoreArray(xx,&x)); 111 PetscCall(PetscLogFlops(2.0*(a->nz) - A->cmap->n)); 112 PetscFunctionReturn(0); 113 } 114