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