1*2c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h> 2*2c733ed4SBarry Smith 3*2c733ed4SBarry Smith PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx) 4*2c733ed4SBarry Smith { 5*2c733ed4SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 6*2c733ed4SBarry Smith PetscErrorCode ierr; 7*2c733ed4SBarry Smith const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j; 8*2c733ed4SBarry Smith PetscInt i,nz,idx,idt,oidx; 9*2c733ed4SBarry Smith const MatScalar *aa=a->a,*v; 10*2c733ed4SBarry Smith PetscScalar s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x; 11*2c733ed4SBarry Smith 12*2c733ed4SBarry Smith PetscFunctionBegin; 13*2c733ed4SBarry Smith ierr = VecCopy(bb,xx);CHKERRQ(ierr); 14*2c733ed4SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 15*2c733ed4SBarry Smith 16*2c733ed4SBarry Smith /* forward solve the U^T */ 17*2c733ed4SBarry Smith idx = 0; 18*2c733ed4SBarry Smith for (i=0; i<n; i++) { 19*2c733ed4SBarry Smith 20*2c733ed4SBarry Smith v = aa + 25*diag[i]; 21*2c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 22*2c733ed4SBarry Smith x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx]; 23*2c733ed4SBarry Smith s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 24*2c733ed4SBarry Smith s2 = v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 25*2c733ed4SBarry Smith s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 26*2c733ed4SBarry Smith s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 27*2c733ed4SBarry Smith s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 28*2c733ed4SBarry Smith v += 25; 29*2c733ed4SBarry Smith 30*2c733ed4SBarry Smith vi = aj + diag[i] + 1; 31*2c733ed4SBarry Smith nz = ai[i+1] - diag[i] - 1; 32*2c733ed4SBarry Smith while (nz--) { 33*2c733ed4SBarry Smith oidx = 5*(*vi++); 34*2c733ed4SBarry Smith x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5; 35*2c733ed4SBarry Smith x[oidx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5; 36*2c733ed4SBarry Smith x[oidx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5; 37*2c733ed4SBarry Smith x[oidx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5; 38*2c733ed4SBarry Smith x[oidx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5; 39*2c733ed4SBarry Smith v += 25; 40*2c733ed4SBarry Smith } 41*2c733ed4SBarry Smith x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;x[3+idx] = s4; x[4+idx] = s5; 42*2c733ed4SBarry Smith idx += 5; 43*2c733ed4SBarry Smith } 44*2c733ed4SBarry Smith /* backward solve the L^T */ 45*2c733ed4SBarry Smith for (i=n-1; i>=0; i--) { 46*2c733ed4SBarry Smith v = aa + 25*diag[i] - 25; 47*2c733ed4SBarry Smith vi = aj + diag[i] - 1; 48*2c733ed4SBarry Smith nz = diag[i] - ai[i]; 49*2c733ed4SBarry Smith idt = 5*i; 50*2c733ed4SBarry Smith s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt];s4 = x[3+idt]; s5 = x[4+idt]; 51*2c733ed4SBarry Smith while (nz--) { 52*2c733ed4SBarry Smith idx = 5*(*vi--); 53*2c733ed4SBarry Smith x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5; 54*2c733ed4SBarry Smith x[idx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5; 55*2c733ed4SBarry Smith x[idx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5; 56*2c733ed4SBarry Smith x[idx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5; 57*2c733ed4SBarry Smith x[idx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5; 58*2c733ed4SBarry Smith v -= 25; 59*2c733ed4SBarry Smith } 60*2c733ed4SBarry Smith } 61*2c733ed4SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 62*2c733ed4SBarry Smith ierr = PetscLogFlops(2.0*25*(a->nz) - 5.0*A->cmap->n);CHKERRQ(ierr); 63*2c733ed4SBarry Smith PetscFunctionReturn(0); 64*2c733ed4SBarry Smith } 65*2c733ed4SBarry Smith 66*2c733ed4SBarry Smith PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat A,Vec bb,Vec xx) 67*2c733ed4SBarry Smith { 68*2c733ed4SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 69*2c733ed4SBarry Smith PetscErrorCode ierr; 70*2c733ed4SBarry Smith const PetscInt n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag; 71*2c733ed4SBarry Smith PetscInt nz,idx,idt,j,i,oidx; 72*2c733ed4SBarry Smith const PetscInt bs =A->rmap->bs,bs2=a->bs2; 73*2c733ed4SBarry Smith const MatScalar *aa=a->a,*v; 74*2c733ed4SBarry Smith PetscScalar s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x; 75*2c733ed4SBarry Smith 76*2c733ed4SBarry Smith PetscFunctionBegin; 77*2c733ed4SBarry Smith ierr = VecCopy(bb,xx);CHKERRQ(ierr); 78*2c733ed4SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 79*2c733ed4SBarry Smith 80*2c733ed4SBarry Smith /* forward solve the U^T */ 81*2c733ed4SBarry Smith idx = 0; 82*2c733ed4SBarry Smith for (i=0; i<n; i++) { 83*2c733ed4SBarry Smith v = aa + bs2*diag[i]; 84*2c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 85*2c733ed4SBarry Smith x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; 86*2c733ed4SBarry Smith x5 = x[4+idx]; 87*2c733ed4SBarry Smith s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 88*2c733ed4SBarry Smith s2 = v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 89*2c733ed4SBarry Smith s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 90*2c733ed4SBarry Smith s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 91*2c733ed4SBarry Smith s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 92*2c733ed4SBarry Smith v -= bs2; 93*2c733ed4SBarry Smith 94*2c733ed4SBarry Smith vi = aj + diag[i] - 1; 95*2c733ed4SBarry Smith nz = diag[i] - diag[i+1] - 1; 96*2c733ed4SBarry Smith for (j=0; j>-nz; j--) { 97*2c733ed4SBarry Smith oidx = bs*vi[j]; 98*2c733ed4SBarry Smith x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5; 99*2c733ed4SBarry Smith x[oidx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5; 100*2c733ed4SBarry Smith x[oidx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5; 101*2c733ed4SBarry Smith x[oidx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5; 102*2c733ed4SBarry Smith x[oidx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5; 103*2c733ed4SBarry Smith v -= bs2; 104*2c733ed4SBarry Smith } 105*2c733ed4SBarry Smith x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3; x[3+idx] = s4; x[4+idx] = s5; 106*2c733ed4SBarry Smith idx += bs; 107*2c733ed4SBarry Smith } 108*2c733ed4SBarry Smith /* backward solve the L^T */ 109*2c733ed4SBarry Smith for (i=n-1; i>=0; i--) { 110*2c733ed4SBarry Smith v = aa + bs2*ai[i]; 111*2c733ed4SBarry Smith vi = aj + ai[i]; 112*2c733ed4SBarry Smith nz = ai[i+1] - ai[i]; 113*2c733ed4SBarry Smith idt = bs*i; 114*2c733ed4SBarry Smith s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt]; s4 = x[3+idt]; s5 = x[4+idt]; 115*2c733ed4SBarry Smith for (j=0; j<nz; j++) { 116*2c733ed4SBarry Smith idx = bs*vi[j]; 117*2c733ed4SBarry Smith x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5; 118*2c733ed4SBarry Smith x[idx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5; 119*2c733ed4SBarry Smith x[idx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5; 120*2c733ed4SBarry Smith x[idx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5; 121*2c733ed4SBarry Smith x[idx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5; 122*2c733ed4SBarry Smith v += bs2; 123*2c733ed4SBarry Smith } 124*2c733ed4SBarry Smith } 125*2c733ed4SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 126*2c733ed4SBarry Smith ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr); 127*2c733ed4SBarry Smith PetscFunctionReturn(0); 128*2c733ed4SBarry Smith } 129