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