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