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 MatSolve_SeqBAIJ_7_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx) 5*2c733ed4SBarry Smith { 6*2c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 7*2c733ed4SBarry Smith const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j; 8*2c733ed4SBarry Smith PetscErrorCode ierr; 9*2c733ed4SBarry Smith PetscInt i,nz,idx,idt,jdx; 10*2c733ed4SBarry Smith const MatScalar *aa=a->a,*v; 11*2c733ed4SBarry Smith PetscScalar *x,s1,s2,s3,s4,s5,s6,s7,x1,x2,x3,x4,x5,x6,x7; 12*2c733ed4SBarry Smith const PetscScalar *b; 13*2c733ed4SBarry Smith 14*2c733ed4SBarry Smith PetscFunctionBegin; 15*2c733ed4SBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 16*2c733ed4SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 17*2c733ed4SBarry Smith /* forward solve the lower triangular */ 18*2c733ed4SBarry Smith idx = 0; 19*2c733ed4SBarry Smith x[0] = b[idx]; x[1] = b[1+idx]; x[2] = b[2+idx]; 20*2c733ed4SBarry Smith x[3] = b[3+idx]; x[4] = b[4+idx]; x[5] = b[5+idx]; 21*2c733ed4SBarry Smith x[6] = b[6+idx]; 22*2c733ed4SBarry Smith for (i=1; i<n; i++) { 23*2c733ed4SBarry Smith v = aa + 49*ai[i]; 24*2c733ed4SBarry Smith vi = aj + ai[i]; 25*2c733ed4SBarry Smith nz = diag[i] - ai[i]; 26*2c733ed4SBarry Smith idx = 7*i; 27*2c733ed4SBarry Smith s1 = b[idx]; s2 = b[1+idx]; s3 = b[2+idx]; 28*2c733ed4SBarry Smith s4 = b[3+idx]; s5 = b[4+idx]; s6 = b[5+idx]; 29*2c733ed4SBarry Smith s7 = b[6+idx]; 30*2c733ed4SBarry Smith while (nz--) { 31*2c733ed4SBarry Smith jdx = 7*(*vi++); 32*2c733ed4SBarry Smith x1 = x[jdx]; x2 = x[1+jdx]; x3 = x[2+jdx]; 33*2c733ed4SBarry Smith x4 = x[3+jdx]; x5 = x[4+jdx]; x6 = x[5+jdx]; 34*2c733ed4SBarry Smith x7 = x[6+jdx]; 35*2c733ed4SBarry Smith s1 -= v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 36*2c733ed4SBarry Smith s2 -= v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 37*2c733ed4SBarry Smith s3 -= v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 38*2c733ed4SBarry Smith s4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 39*2c733ed4SBarry Smith s5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 40*2c733ed4SBarry Smith s6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 41*2c733ed4SBarry Smith s7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 42*2c733ed4SBarry Smith v += 49; 43*2c733ed4SBarry Smith } 44*2c733ed4SBarry Smith x[idx] = s1; 45*2c733ed4SBarry Smith x[1+idx] = s2; 46*2c733ed4SBarry Smith x[2+idx] = s3; 47*2c733ed4SBarry Smith x[3+idx] = s4; 48*2c733ed4SBarry Smith x[4+idx] = s5; 49*2c733ed4SBarry Smith x[5+idx] = s6; 50*2c733ed4SBarry Smith x[6+idx] = s7; 51*2c733ed4SBarry Smith } 52*2c733ed4SBarry Smith /* backward solve the upper triangular */ 53*2c733ed4SBarry Smith for (i=n-1; i>=0; i--) { 54*2c733ed4SBarry Smith v = aa + 49*diag[i] + 49; 55*2c733ed4SBarry Smith vi = aj + diag[i] + 1; 56*2c733ed4SBarry Smith nz = ai[i+1] - diag[i] - 1; 57*2c733ed4SBarry Smith idt = 7*i; 58*2c733ed4SBarry Smith s1 = x[idt]; s2 = x[1+idt]; 59*2c733ed4SBarry Smith s3 = x[2+idt]; s4 = x[3+idt]; 60*2c733ed4SBarry Smith s5 = x[4+idt]; s6 = x[5+idt]; 61*2c733ed4SBarry Smith s7 = x[6+idt]; 62*2c733ed4SBarry Smith while (nz--) { 63*2c733ed4SBarry Smith idx = 7*(*vi++); 64*2c733ed4SBarry Smith x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; 65*2c733ed4SBarry Smith x4 = x[3+idx]; x5 = x[4+idx]; x6 = x[5+idx]; 66*2c733ed4SBarry Smith x7 = x[6+idx]; 67*2c733ed4SBarry Smith s1 -= v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 68*2c733ed4SBarry Smith s2 -= v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 69*2c733ed4SBarry Smith s3 -= v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 70*2c733ed4SBarry Smith s4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 71*2c733ed4SBarry Smith s5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 72*2c733ed4SBarry Smith s6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 73*2c733ed4SBarry Smith s7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 74*2c733ed4SBarry Smith v += 49; 75*2c733ed4SBarry Smith } 76*2c733ed4SBarry Smith v = aa + 49*diag[i]; 77*2c733ed4SBarry Smith x[idt] = v[0]*s1 + v[7]*s2 + v[14]*s3 + v[21]*s4 78*2c733ed4SBarry Smith + v[28]*s5 + v[35]*s6 + v[42]*s7; 79*2c733ed4SBarry Smith x[1+idt] = v[1]*s1 + v[8]*s2 + v[15]*s3 + v[22]*s4 80*2c733ed4SBarry Smith + v[29]*s5 + v[36]*s6 + v[43]*s7; 81*2c733ed4SBarry Smith x[2+idt] = v[2]*s1 + v[9]*s2 + v[16]*s3 + v[23]*s4 82*2c733ed4SBarry Smith + v[30]*s5 + v[37]*s6 + v[44]*s7; 83*2c733ed4SBarry Smith x[3+idt] = v[3]*s1 + v[10]*s2 + v[17]*s3 + v[24]*s4 84*2c733ed4SBarry Smith + v[31]*s5 + v[38]*s6 + v[45]*s7; 85*2c733ed4SBarry Smith x[4+idt] = v[4]*s1 + v[11]*s2 + v[18]*s3 + v[25]*s4 86*2c733ed4SBarry Smith + v[32]*s5 + v[39]*s6 + v[46]*s7; 87*2c733ed4SBarry Smith x[5+idt] = v[5]*s1 + v[12]*s2 + v[19]*s3 + v[26]*s4 88*2c733ed4SBarry Smith + v[33]*s5 + v[40]*s6 + v[47]*s7; 89*2c733ed4SBarry Smith x[6+idt] = v[6]*s1 + v[13]*s2 + v[20]*s3 + v[27]*s4 90*2c733ed4SBarry Smith + v[34]*s5 + v[41]*s6 + v[48]*s7; 91*2c733ed4SBarry Smith } 92*2c733ed4SBarry Smith 93*2c733ed4SBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 94*2c733ed4SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 95*2c733ed4SBarry Smith ierr = PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n);CHKERRQ(ierr); 96*2c733ed4SBarry Smith PetscFunctionReturn(0); 97*2c733ed4SBarry Smith } 98*2c733ed4SBarry Smith 99*2c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering(Mat A,Vec bb,Vec xx) 100*2c733ed4SBarry Smith { 101*2c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 102*2c733ed4SBarry Smith const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*adiag=a->diag; 103*2c733ed4SBarry Smith PetscErrorCode ierr; 104*2c733ed4SBarry Smith PetscInt i,k,nz,idx,jdx,idt; 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 *x; 108*2c733ed4SBarry Smith const PetscScalar *b; 109*2c733ed4SBarry Smith PetscScalar s1,s2,s3,s4,s5,s6,s7,x1,x2,x3,x4,x5,x6,x7; 110*2c733ed4SBarry Smith 111*2c733ed4SBarry Smith PetscFunctionBegin; 112*2c733ed4SBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 113*2c733ed4SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 114*2c733ed4SBarry Smith /* forward solve the lower triangular */ 115*2c733ed4SBarry Smith idx = 0; 116*2c733ed4SBarry Smith x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx];x[3] = b[3+idx]; 117*2c733ed4SBarry Smith x[4] = b[4+idx];x[5] = b[5+idx];x[6] = b[6+idx]; 118*2c733ed4SBarry Smith for (i=1; i<n; i++) { 119*2c733ed4SBarry Smith v = aa + bs2*ai[i]; 120*2c733ed4SBarry Smith vi = aj + ai[i]; 121*2c733ed4SBarry Smith nz = ai[i+1] - ai[i]; 122*2c733ed4SBarry Smith idx = bs*i; 123*2c733ed4SBarry Smith s1 = b[idx];s2 = b[1+idx];s3 = b[2+idx];s4 = b[3+idx]; 124*2c733ed4SBarry Smith s5 = b[4+idx];s6 = b[5+idx];s7 = b[6+idx]; 125*2c733ed4SBarry Smith for (k=0; k<nz; k++) { 126*2c733ed4SBarry Smith jdx = bs*vi[k]; 127*2c733ed4SBarry Smith x1 = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx];x4 =x[3+jdx]; 128*2c733ed4SBarry Smith x5 = x[4+jdx]; x6 = x[5+jdx];x7 = x[6+jdx]; 129*2c733ed4SBarry Smith s1 -= v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 130*2c733ed4SBarry Smith s2 -= v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 131*2c733ed4SBarry Smith s3 -= v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 132*2c733ed4SBarry Smith s4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 133*2c733ed4SBarry Smith s5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 134*2c733ed4SBarry Smith s6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 135*2c733ed4SBarry Smith s7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 136*2c733ed4SBarry Smith v += bs2; 137*2c733ed4SBarry Smith } 138*2c733ed4SBarry Smith 139*2c733ed4SBarry Smith x[idx] = s1; 140*2c733ed4SBarry Smith x[1+idx] = s2; 141*2c733ed4SBarry Smith x[2+idx] = s3; 142*2c733ed4SBarry Smith x[3+idx] = s4; 143*2c733ed4SBarry Smith x[4+idx] = s5; 144*2c733ed4SBarry Smith x[5+idx] = s6; 145*2c733ed4SBarry Smith x[6+idx] = s7; 146*2c733ed4SBarry Smith } 147*2c733ed4SBarry Smith 148*2c733ed4SBarry Smith /* backward solve the upper triangular */ 149*2c733ed4SBarry Smith for (i=n-1; i>=0; i--) { 150*2c733ed4SBarry Smith v = aa + bs2*(adiag[i+1]+1); 151*2c733ed4SBarry Smith vi = aj + adiag[i+1]+1; 152*2c733ed4SBarry Smith nz = adiag[i] - adiag[i+1]-1; 153*2c733ed4SBarry Smith idt = bs*i; 154*2c733ed4SBarry Smith s1 = x[idt]; s2 = x[1+idt];s3 = x[2+idt];s4 = x[3+idt]; 155*2c733ed4SBarry Smith s5 = x[4+idt];s6 = x[5+idt];s7 = x[6+idt]; 156*2c733ed4SBarry Smith for (k=0; k<nz; k++) { 157*2c733ed4SBarry Smith idx = bs*vi[k]; 158*2c733ed4SBarry Smith x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx];x4 = x[3+idx]; 159*2c733ed4SBarry Smith x5 = x[4+idx];x6 = x[5+idx];x7 = x[6+idx]; 160*2c733ed4SBarry Smith s1 -= v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 161*2c733ed4SBarry Smith s2 -= v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 162*2c733ed4SBarry Smith s3 -= v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 163*2c733ed4SBarry Smith s4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 164*2c733ed4SBarry Smith s5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 165*2c733ed4SBarry Smith s6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 166*2c733ed4SBarry Smith s7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 167*2c733ed4SBarry Smith v += bs2; 168*2c733ed4SBarry Smith } 169*2c733ed4SBarry Smith /* x = inv_diagonal*x */ 170*2c733ed4SBarry Smith x[idt] = v[0]*s1 + v[7]*s2 + v[14]*s3 + v[21]*s4 + v[28]*s5 + v[35]*s6 + v[42]*s7; 171*2c733ed4SBarry Smith x[1+idt] = v[1]*s1 + v[8]*s2 + v[15]*s3 + v[22]*s4 + v[29]*s5 + v[36]*s6 + v[43]*s7; 172*2c733ed4SBarry Smith x[2+idt] = v[2]*s1 + v[9]*s2 + v[16]*s3 + v[23]*s4 + v[30]*s5 + v[37]*s6 + v[44]*s7; 173*2c733ed4SBarry Smith x[3+idt] = v[3]*s1 + v[10]*s2 + v[17]*s3 + v[24]*s4 + v[31]*s5 + v[38]*s6 + v[45]*s7; 174*2c733ed4SBarry Smith x[4+idt] = v[4]*s1 + v[11]*s2 + v[18]*s3 + v[25]*s4 + v[32]*s5 + v[39]*s6 + v[46]*s7; 175*2c733ed4SBarry Smith x[5+idt] = v[5]*s1 + v[12]*s2 + v[19]*s3 + v[26]*s4 + v[33]*s5 + v[40]*s6 + v[47]*s7; 176*2c733ed4SBarry Smith x[6+idt] = v[6]*s1 + v[13]*s2 + v[20]*s3 + v[27]*s4 + v[34]*s5 + v[41]*s6 + v[48]*s7; 177*2c733ed4SBarry Smith } 178*2c733ed4SBarry Smith 179*2c733ed4SBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 180*2c733ed4SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 181*2c733ed4SBarry Smith ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr); 182*2c733ed4SBarry Smith PetscFunctionReturn(0); 183*2c733ed4SBarry Smith } 184*2c733ed4SBarry Smith 185