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_4_inplace(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 *r,*c,*rout,*cout; 10*2c733ed4SBarry Smith const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j; 11*2c733ed4SBarry Smith PetscInt i,nz,idx,idt,ii,ic,ir,oidx; 12*2c733ed4SBarry Smith const MatScalar *aa=a->a,*v; 13*2c733ed4SBarry Smith PetscScalar s1,s2,s3,s4,x1,x2,x3,x4,*x,*t; 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 t = 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 ii = 0; 26*2c733ed4SBarry Smith for (i=0; i<n; i++) { 27*2c733ed4SBarry Smith ic = 4*c[i]; 28*2c733ed4SBarry Smith t[ii] = b[ic]; 29*2c733ed4SBarry Smith t[ii+1] = b[ic+1]; 30*2c733ed4SBarry Smith t[ii+2] = b[ic+2]; 31*2c733ed4SBarry Smith t[ii+3] = b[ic+3]; 32*2c733ed4SBarry Smith ii += 4; 33*2c733ed4SBarry Smith } 34*2c733ed4SBarry Smith 35*2c733ed4SBarry Smith /* forward solve the U^T */ 36*2c733ed4SBarry Smith idx = 0; 37*2c733ed4SBarry Smith for (i=0; i<n; i++) { 38*2c733ed4SBarry Smith 39*2c733ed4SBarry Smith v = aa + 16*diag[i]; 40*2c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 41*2c733ed4SBarry Smith x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx]; 42*2c733ed4SBarry Smith s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 43*2c733ed4SBarry Smith s2 = v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 44*2c733ed4SBarry Smith s3 = v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 45*2c733ed4SBarry Smith s4 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 46*2c733ed4SBarry Smith v += 16; 47*2c733ed4SBarry Smith 48*2c733ed4SBarry Smith vi = aj + diag[i] + 1; 49*2c733ed4SBarry Smith nz = ai[i+1] - diag[i] - 1; 50*2c733ed4SBarry Smith while (nz--) { 51*2c733ed4SBarry Smith oidx = 4*(*vi++); 52*2c733ed4SBarry Smith t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4; 53*2c733ed4SBarry Smith t[oidx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4; 54*2c733ed4SBarry Smith t[oidx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4; 55*2c733ed4SBarry Smith t[oidx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4; 56*2c733ed4SBarry Smith v += 16; 57*2c733ed4SBarry Smith } 58*2c733ed4SBarry Smith t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4; 59*2c733ed4SBarry Smith idx += 4; 60*2c733ed4SBarry Smith } 61*2c733ed4SBarry Smith /* backward solve the L^T */ 62*2c733ed4SBarry Smith for (i=n-1; i>=0; i--) { 63*2c733ed4SBarry Smith v = aa + 16*diag[i] - 16; 64*2c733ed4SBarry Smith vi = aj + diag[i] - 1; 65*2c733ed4SBarry Smith nz = diag[i] - ai[i]; 66*2c733ed4SBarry Smith idt = 4*i; 67*2c733ed4SBarry Smith s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt]; 68*2c733ed4SBarry Smith while (nz--) { 69*2c733ed4SBarry Smith idx = 4*(*vi--); 70*2c733ed4SBarry Smith t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4; 71*2c733ed4SBarry Smith t[idx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4; 72*2c733ed4SBarry Smith t[idx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4; 73*2c733ed4SBarry Smith t[idx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4; 74*2c733ed4SBarry Smith v -= 16; 75*2c733ed4SBarry Smith } 76*2c733ed4SBarry Smith } 77*2c733ed4SBarry Smith 78*2c733ed4SBarry Smith /* copy t into x according to permutation */ 79*2c733ed4SBarry Smith ii = 0; 80*2c733ed4SBarry Smith for (i=0; i<n; i++) { 81*2c733ed4SBarry Smith ir = 4*r[i]; 82*2c733ed4SBarry Smith x[ir] = t[ii]; 83*2c733ed4SBarry Smith x[ir+1] = t[ii+1]; 84*2c733ed4SBarry Smith x[ir+2] = t[ii+2]; 85*2c733ed4SBarry Smith x[ir+3] = t[ii+3]; 86*2c733ed4SBarry Smith ii += 4; 87*2c733ed4SBarry Smith } 88*2c733ed4SBarry Smith 89*2c733ed4SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 90*2c733ed4SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 91*2c733ed4SBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 92*2c733ed4SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 93*2c733ed4SBarry Smith ierr = PetscLogFlops(2.0*16*(a->nz) - 4.0*A->cmap->n);CHKERRQ(ierr); 94*2c733ed4SBarry Smith PetscFunctionReturn(0); 95*2c733ed4SBarry Smith } 96*2c733ed4SBarry Smith 97*2c733ed4SBarry Smith PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat A,Vec bb,Vec xx) 98*2c733ed4SBarry Smith { 99*2c733ed4SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 100*2c733ed4SBarry Smith PetscErrorCode ierr; 101*2c733ed4SBarry Smith IS iscol=a->col,isrow=a->row; 102*2c733ed4SBarry Smith const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag; 103*2c733ed4SBarry Smith const PetscInt *r,*c,*rout,*cout; 104*2c733ed4SBarry Smith PetscInt nz,idx,idt,j,i,oidx,ii,ic,ir; 105*2c733ed4SBarry Smith const PetscInt bs =A->rmap->bs,bs2=a->bs2; 106*2c733ed4SBarry Smith const MatScalar *aa=a->a,*v; 107*2c733ed4SBarry Smith PetscScalar s1,s2,s3,s4,x1,x2,x3,x4,*x,*t; 108*2c733ed4SBarry Smith const PetscScalar *b; 109*2c733ed4SBarry Smith 110*2c733ed4SBarry Smith PetscFunctionBegin; 111*2c733ed4SBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 112*2c733ed4SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 113*2c733ed4SBarry Smith t = a->solve_work; 114*2c733ed4SBarry Smith 115*2c733ed4SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 116*2c733ed4SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 117*2c733ed4SBarry Smith 118*2c733ed4SBarry Smith /* copy b into temp work space according to permutation */ 119*2c733ed4SBarry Smith for (i=0; i<n; i++) { 120*2c733ed4SBarry Smith ii = bs*i; ic = bs*c[i]; 121*2c733ed4SBarry Smith t[ii] = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3]; 122*2c733ed4SBarry Smith } 123*2c733ed4SBarry Smith 124*2c733ed4SBarry Smith /* forward solve the U^T */ 125*2c733ed4SBarry Smith idx = 0; 126*2c733ed4SBarry Smith for (i=0; i<n; i++) { 127*2c733ed4SBarry Smith v = aa + bs2*diag[i]; 128*2c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 129*2c733ed4SBarry Smith x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx]; 130*2c733ed4SBarry Smith s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 131*2c733ed4SBarry Smith s2 = v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 132*2c733ed4SBarry Smith s3 = v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 133*2c733ed4SBarry Smith s4 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 134*2c733ed4SBarry Smith v -= bs2; 135*2c733ed4SBarry Smith 136*2c733ed4SBarry Smith vi = aj + diag[i] - 1; 137*2c733ed4SBarry Smith nz = diag[i] - diag[i+1] - 1; 138*2c733ed4SBarry Smith for (j=0; j>-nz; j--) { 139*2c733ed4SBarry Smith oidx = bs*vi[j]; 140*2c733ed4SBarry Smith t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4; 141*2c733ed4SBarry Smith t[oidx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4; 142*2c733ed4SBarry Smith t[oidx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4; 143*2c733ed4SBarry Smith t[oidx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4; 144*2c733ed4SBarry Smith v -= bs2; 145*2c733ed4SBarry Smith } 146*2c733ed4SBarry Smith t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3; t[3+idx] = s4; 147*2c733ed4SBarry Smith idx += bs; 148*2c733ed4SBarry Smith } 149*2c733ed4SBarry Smith /* backward solve the L^T */ 150*2c733ed4SBarry Smith for (i=n-1; i>=0; i--) { 151*2c733ed4SBarry Smith v = aa + bs2*ai[i]; 152*2c733ed4SBarry Smith vi = aj + ai[i]; 153*2c733ed4SBarry Smith nz = ai[i+1] - ai[i]; 154*2c733ed4SBarry Smith idt = bs*i; 155*2c733ed4SBarry Smith s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt]; s4 = t[3+idt]; 156*2c733ed4SBarry Smith for (j=0; j<nz; j++) { 157*2c733ed4SBarry Smith idx = bs*vi[j]; 158*2c733ed4SBarry Smith t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4; 159*2c733ed4SBarry Smith t[idx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4; 160*2c733ed4SBarry Smith t[idx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4; 161*2c733ed4SBarry Smith t[idx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4; 162*2c733ed4SBarry Smith v += bs2; 163*2c733ed4SBarry Smith } 164*2c733ed4SBarry Smith } 165*2c733ed4SBarry Smith 166*2c733ed4SBarry Smith /* copy t into x according to permutation */ 167*2c733ed4SBarry Smith for (i=0; i<n; i++) { 168*2c733ed4SBarry Smith ii = bs*i; ir = bs*r[i]; 169*2c733ed4SBarry Smith x[ir] = t[ii]; x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2]; x[ir+3] = t[ii+3]; 170*2c733ed4SBarry Smith } 171*2c733ed4SBarry Smith 172*2c733ed4SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 173*2c733ed4SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 174*2c733ed4SBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 175*2c733ed4SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 176*2c733ed4SBarry Smith ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr); 177*2c733ed4SBarry Smith PetscFunctionReturn(0); 178*2c733ed4SBarry Smith } 179