1*4e2b4712SSatish Balay #ifdef PETSC_RCS_HEADER 2*4e2b4712SSatish Balay static char vcid[] = "$Id: baijfact.c,v 1.53 1997/10/19 03:26:08 bsmith Exp balay $"; 3*4e2b4712SSatish Balay #endif 4*4e2b4712SSatish Balay /* 5*4e2b4712SSatish Balay Factorization code for BAIJ format. 6*4e2b4712SSatish Balay */ 7*4e2b4712SSatish Balay 8*4e2b4712SSatish Balay #include "src/mat/impls/baij/seq/baij.h" 9*4e2b4712SSatish Balay #include "src/vec/vecimpl.h" 10*4e2b4712SSatish Balay #include "src/inline/ilu.h" 11*4e2b4712SSatish Balay 12*4e2b4712SSatish Balay /* ----------------------------------------------------------- */ 13*4e2b4712SSatish Balay #undef __FUNC__ 14*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_N" 15*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx) 16*4e2b4712SSatish Balay { 17*4e2b4712SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 18*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 19*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j; 20*4e2b4712SSatish Balay int nz,bs=a->bs,bs2=a->bs2,*rout,*cout; 21*4e2b4712SSatish Balay Scalar *aa=a->a,*sum; 22*4e2b4712SSatish Balay register Scalar *x,*b,*lsum,*tmp,*v; 23*4e2b4712SSatish Balay 24*4e2b4712SSatish Balay PetscFunctionBegin; 25*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 26*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 27*4e2b4712SSatish Balay tmp = a->solve_work; 28*4e2b4712SSatish Balay 29*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 30*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 31*4e2b4712SSatish Balay 32*4e2b4712SSatish Balay /* forward solve the lower triangular */ 33*4e2b4712SSatish Balay PetscMemcpy(tmp,b + bs*(*r++), bs*sizeof(Scalar)); 34*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 35*4e2b4712SSatish Balay v = aa + bs2*ai[i]; 36*4e2b4712SSatish Balay vi = aj + ai[i]; 37*4e2b4712SSatish Balay nz = a->diag[i] - ai[i]; 38*4e2b4712SSatish Balay sum = tmp + bs*i; 39*4e2b4712SSatish Balay PetscMemcpy(sum,b+bs*(*r++),bs*sizeof(Scalar)); 40*4e2b4712SSatish Balay while (nz--) { 41*4e2b4712SSatish Balay Kernel_v_gets_v_minus_A_times_w(bs,sum,v,tmp+bs*(*vi++)); 42*4e2b4712SSatish Balay v += bs2; 43*4e2b4712SSatish Balay } 44*4e2b4712SSatish Balay } 45*4e2b4712SSatish Balay /* backward solve the upper triangular */ 46*4e2b4712SSatish Balay lsum = a->solve_work + a->n; 47*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 48*4e2b4712SSatish Balay v = aa + bs2*(a->diag[i] + 1); 49*4e2b4712SSatish Balay vi = aj + a->diag[i] + 1; 50*4e2b4712SSatish Balay nz = ai[i+1] - a->diag[i] - 1; 51*4e2b4712SSatish Balay PetscMemcpy(lsum,tmp+i*bs,bs*sizeof(Scalar)); 52*4e2b4712SSatish Balay while (nz--) { 53*4e2b4712SSatish Balay Kernel_v_gets_v_minus_A_times_w(bs,lsum,v,tmp+bs*(*vi++)); 54*4e2b4712SSatish Balay v += bs2; 55*4e2b4712SSatish Balay } 56*4e2b4712SSatish Balay Kernel_w_gets_A_times_v(bs,lsum,aa+bs2*a->diag[i],tmp+i*bs); 57*4e2b4712SSatish Balay PetscMemcpy(x + bs*(*c--),tmp+i*bs,bs*sizeof(Scalar)); 58*4e2b4712SSatish Balay } 59*4e2b4712SSatish Balay 60*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 61*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 62*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 63*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 64*4e2b4712SSatish Balay PLogFlops(2*(a->bs2)*(a->nz) - a->n); 65*4e2b4712SSatish Balay PetscFunctionReturn(0); 66*4e2b4712SSatish Balay } 67*4e2b4712SSatish Balay 68*4e2b4712SSatish Balay #undef __FUNC__ 69*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_7" 70*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_7(Mat A,Vec bb,Vec xx) 71*4e2b4712SSatish Balay { 72*4e2b4712SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 73*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 74*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout; 75*4e2b4712SSatish Balay int *diag = a->diag; 76*4e2b4712SSatish Balay Scalar *aa=a->a,sum1,sum2,sum3,sum4,sum5,sum6,sum7,x1,x2,x3,x4,x5,x6,x7; 77*4e2b4712SSatish Balay register Scalar *x,*b,*tmp,*v; 78*4e2b4712SSatish Balay 79*4e2b4712SSatish Balay PetscFunctionBegin; 80*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 81*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 82*4e2b4712SSatish Balay tmp = a->solve_work; 83*4e2b4712SSatish Balay 84*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 85*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 86*4e2b4712SSatish Balay 87*4e2b4712SSatish Balay /* forward solve the lower triangular */ 88*4e2b4712SSatish Balay idx = 7*(*r++); 89*4e2b4712SSatish Balay tmp[0] = b[idx]; tmp[1] = b[1+idx]; 90*4e2b4712SSatish Balay tmp[2] = b[2+idx]; tmp[3] = b[3+idx]; tmp[4] = b[4+idx]; 91*4e2b4712SSatish Balay tmp[5] = b[5+idx]; tmp[6] = b[6+idx]; 92*4e2b4712SSatish Balay 93*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 94*4e2b4712SSatish Balay v = aa + 49*ai[i]; 95*4e2b4712SSatish Balay vi = aj + ai[i]; 96*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 97*4e2b4712SSatish Balay idx = 7*(*r++); 98*4e2b4712SSatish Balay sum1 = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx]; 99*4e2b4712SSatish Balay sum5 = b[4+idx];sum6 = b[5+idx];sum7 = b[6+idx]; 100*4e2b4712SSatish Balay while (nz--) { 101*4e2b4712SSatish Balay idx = 7*(*vi++); 102*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx];x3 = tmp[2+idx]; 103*4e2b4712SSatish Balay x4 = tmp[3+idx];x5 = tmp[4+idx]; 104*4e2b4712SSatish Balay x6 = tmp[5+idx];x7 = tmp[6+idx]; 105*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 106*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 107*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 108*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 109*4e2b4712SSatish Balay sum5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 110*4e2b4712SSatish Balay sum6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 111*4e2b4712SSatish Balay sum7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 112*4e2b4712SSatish Balay v += 49; 113*4e2b4712SSatish Balay } 114*4e2b4712SSatish Balay idx = 7*i; 115*4e2b4712SSatish Balay tmp[idx] = sum1;tmp[1+idx] = sum2; 116*4e2b4712SSatish Balay tmp[2+idx] = sum3;tmp[3+idx] = sum4; tmp[4+idx] = sum5; 117*4e2b4712SSatish Balay tmp[5+idx] = sum6;tmp[6+idx] = sum7; 118*4e2b4712SSatish Balay } 119*4e2b4712SSatish Balay /* backward solve the upper triangular */ 120*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 121*4e2b4712SSatish Balay v = aa + 49*diag[i] + 49; 122*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 123*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 124*4e2b4712SSatish Balay idt = 7*i; 125*4e2b4712SSatish Balay sum1 = tmp[idt]; sum2 = tmp[1+idt]; 126*4e2b4712SSatish Balay sum3 = tmp[2+idt];sum4 = tmp[3+idt]; sum5 = tmp[4+idt]; 127*4e2b4712SSatish Balay sum6 = tmp[5+idt];sum7 = tmp[6+idt]; 128*4e2b4712SSatish Balay while (nz--) { 129*4e2b4712SSatish Balay idx = 7*(*vi++); 130*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; 131*4e2b4712SSatish Balay x3 = tmp[2+idx]; x4 = tmp[3+idx]; x5 = tmp[4+idx]; 132*4e2b4712SSatish Balay x6 = tmp[5+idx]; x7 = tmp[6+idx]; 133*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 134*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 135*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 136*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 137*4e2b4712SSatish Balay sum5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 138*4e2b4712SSatish Balay sum6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 139*4e2b4712SSatish Balay sum7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 140*4e2b4712SSatish Balay v += 49; 141*4e2b4712SSatish Balay } 142*4e2b4712SSatish Balay idc = 7*(*c--); 143*4e2b4712SSatish Balay v = aa + 49*diag[i]; 144*4e2b4712SSatish Balay x[idc] = tmp[idt] = v[0]*sum1+v[7]*sum2+v[14]*sum3+ 145*4e2b4712SSatish Balay v[21]*sum4+v[28]*sum5+v[35]*sum6+v[42]*sum7; 146*4e2b4712SSatish Balay x[1+idc] = tmp[1+idt] = v[1]*sum1+v[8]*sum2+v[15]*sum3+ 147*4e2b4712SSatish Balay v[22]*sum4+v[29]*sum5+v[36]*sum6+v[43]*sum7; 148*4e2b4712SSatish Balay x[2+idc] = tmp[2+idt] = v[2]*sum1+v[9]*sum2+v[16]*sum3+ 149*4e2b4712SSatish Balay v[23]*sum4+v[30]*sum5+v[37]*sum6+v[44]*sum7; 150*4e2b4712SSatish Balay x[3+idc] = tmp[3+idt] = v[3]*sum1+v[10]*sum2+v[17]*sum3+ 151*4e2b4712SSatish Balay v[24]*sum4+v[31]*sum5+v[38]*sum6+v[45]*sum7; 152*4e2b4712SSatish Balay x[4+idc] = tmp[4+idt] = v[4]*sum1+v[11]*sum2+v[18]*sum3+ 153*4e2b4712SSatish Balay v[25]*sum4+v[32]*sum5+v[39]*sum6+v[46]*sum7; 154*4e2b4712SSatish Balay x[5+idc] = tmp[5+idt] = v[5]*sum1+v[12]*sum2+v[19]*sum3+ 155*4e2b4712SSatish Balay v[26]*sum4+v[33]*sum5+v[40]*sum6+v[47]*sum7; 156*4e2b4712SSatish Balay x[6+idc] = tmp[6+idt] = v[6]*sum1+v[13]*sum2+v[20]*sum3+ 157*4e2b4712SSatish Balay v[27]*sum4+v[34]*sum5+v[41]*sum6+v[48]*sum7; 158*4e2b4712SSatish Balay } 159*4e2b4712SSatish Balay 160*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 161*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 162*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 163*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 164*4e2b4712SSatish Balay PLogFlops(2*49*(a->nz) - a->n); 165*4e2b4712SSatish Balay PetscFunctionReturn(0); 166*4e2b4712SSatish Balay } 167*4e2b4712SSatish Balay 168*4e2b4712SSatish Balay #undef __FUNC__ 169*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_5" 170*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_5(Mat A,Vec bb,Vec xx) 171*4e2b4712SSatish Balay { 172*4e2b4712SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 173*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 174*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout; 175*4e2b4712SSatish Balay int *diag = a->diag; 176*4e2b4712SSatish Balay Scalar *aa=a->a,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 177*4e2b4712SSatish Balay register Scalar *x,*b,*tmp,*v; 178*4e2b4712SSatish Balay 179*4e2b4712SSatish Balay PetscFunctionBegin; 180*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 181*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 182*4e2b4712SSatish Balay tmp = a->solve_work; 183*4e2b4712SSatish Balay 184*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 185*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 186*4e2b4712SSatish Balay 187*4e2b4712SSatish Balay /* forward solve the lower triangular */ 188*4e2b4712SSatish Balay idx = 5*(*r++); 189*4e2b4712SSatish Balay tmp[0] = b[idx]; tmp[1] = b[1+idx]; 190*4e2b4712SSatish Balay tmp[2] = b[2+idx]; tmp[3] = b[3+idx]; tmp[4] = b[4+idx]; 191*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 192*4e2b4712SSatish Balay v = aa + 25*ai[i]; 193*4e2b4712SSatish Balay vi = aj + ai[i]; 194*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 195*4e2b4712SSatish Balay idx = 5*(*r++); 196*4e2b4712SSatish Balay sum1 = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx]; 197*4e2b4712SSatish Balay sum5 = b[4+idx]; 198*4e2b4712SSatish Balay while (nz--) { 199*4e2b4712SSatish Balay idx = 5*(*vi++); 200*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx];x3 = tmp[2+idx]; 201*4e2b4712SSatish Balay x4 = tmp[3+idx];x5 = tmp[4+idx]; 202*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 203*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 204*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 205*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 206*4e2b4712SSatish Balay sum5 -= v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 207*4e2b4712SSatish Balay v += 25; 208*4e2b4712SSatish Balay } 209*4e2b4712SSatish Balay idx = 5*i; 210*4e2b4712SSatish Balay tmp[idx] = sum1;tmp[1+idx] = sum2; 211*4e2b4712SSatish Balay tmp[2+idx] = sum3;tmp[3+idx] = sum4; tmp[4+idx] = sum5; 212*4e2b4712SSatish Balay } 213*4e2b4712SSatish Balay /* backward solve the upper triangular */ 214*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 215*4e2b4712SSatish Balay v = aa + 25*diag[i] + 25; 216*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 217*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 218*4e2b4712SSatish Balay idt = 5*i; 219*4e2b4712SSatish Balay sum1 = tmp[idt]; sum2 = tmp[1+idt]; 220*4e2b4712SSatish Balay sum3 = tmp[2+idt];sum4 = tmp[3+idt]; sum5 = tmp[4+idt]; 221*4e2b4712SSatish Balay while (nz--) { 222*4e2b4712SSatish Balay idx = 5*(*vi++); 223*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; 224*4e2b4712SSatish Balay x3 = tmp[2+idx]; x4 = tmp[3+idx]; x5 = tmp[4+idx]; 225*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 226*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 227*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 228*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 229*4e2b4712SSatish Balay sum5 -= v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 230*4e2b4712SSatish Balay v += 25; 231*4e2b4712SSatish Balay } 232*4e2b4712SSatish Balay idc = 5*(*c--); 233*4e2b4712SSatish Balay v = aa + 25*diag[i]; 234*4e2b4712SSatish Balay x[idc] = tmp[idt] = v[0]*sum1+v[5]*sum2+v[10]*sum3+ 235*4e2b4712SSatish Balay v[15]*sum4+v[20]*sum5; 236*4e2b4712SSatish Balay x[1+idc] = tmp[1+idt] = v[1]*sum1+v[6]*sum2+v[11]*sum3+ 237*4e2b4712SSatish Balay v[16]*sum4+v[21]*sum5; 238*4e2b4712SSatish Balay x[2+idc] = tmp[2+idt] = v[2]*sum1+v[7]*sum2+v[12]*sum3+ 239*4e2b4712SSatish Balay v[17]*sum4+v[22]*sum5; 240*4e2b4712SSatish Balay x[3+idc] = tmp[3+idt] = v[3]*sum1+v[8]*sum2+v[13]*sum3+ 241*4e2b4712SSatish Balay v[18]*sum4+v[23]*sum5; 242*4e2b4712SSatish Balay x[4+idc] = tmp[4+idt] = v[4]*sum1+v[9]*sum2+v[14]*sum3+ 243*4e2b4712SSatish Balay v[19]*sum4+v[24]*sum5; 244*4e2b4712SSatish Balay } 245*4e2b4712SSatish Balay 246*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 247*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 248*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 249*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 250*4e2b4712SSatish Balay PLogFlops(2*25*(a->nz) - a->n); 251*4e2b4712SSatish Balay PetscFunctionReturn(0); 252*4e2b4712SSatish Balay } 253*4e2b4712SSatish Balay 254*4e2b4712SSatish Balay #undef __FUNC__ 255*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_4" 256*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_4(Mat A,Vec bb,Vec xx) 257*4e2b4712SSatish Balay { 258*4e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 259*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 260*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout; 261*4e2b4712SSatish Balay int *diag = a->diag; 262*4e2b4712SSatish Balay Scalar *aa=a->a,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 263*4e2b4712SSatish Balay register Scalar *x,*b,*tmp,*v; 264*4e2b4712SSatish Balay 265*4e2b4712SSatish Balay PetscFunctionBegin; 266*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 267*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 268*4e2b4712SSatish Balay tmp = a->solve_work; 269*4e2b4712SSatish Balay 270*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 271*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 272*4e2b4712SSatish Balay 273*4e2b4712SSatish Balay /* forward solve the lower triangular */ 274*4e2b4712SSatish Balay idx = 4*(*r++); 275*4e2b4712SSatish Balay tmp[0] = b[idx]; tmp[1] = b[1+idx]; 276*4e2b4712SSatish Balay tmp[2] = b[2+idx]; tmp[3] = b[3+idx]; 277*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 278*4e2b4712SSatish Balay v = aa + 16*ai[i]; 279*4e2b4712SSatish Balay vi = aj + ai[i]; 280*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 281*4e2b4712SSatish Balay idx = 4*(*r++); 282*4e2b4712SSatish Balay sum1 = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx]; 283*4e2b4712SSatish Balay while (nz--) { 284*4e2b4712SSatish Balay idx = 4*(*vi++); 285*4e2b4712SSatish Balay x1 = tmp[idx];x2 = tmp[1+idx];x3 = tmp[2+idx];x4 = tmp[3+idx]; 286*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 287*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 288*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 289*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 290*4e2b4712SSatish Balay v += 16; 291*4e2b4712SSatish Balay } 292*4e2b4712SSatish Balay idx = 4*i; 293*4e2b4712SSatish Balay tmp[idx] = sum1;tmp[1+idx] = sum2; 294*4e2b4712SSatish Balay tmp[2+idx] = sum3;tmp[3+idx] = sum4; 295*4e2b4712SSatish Balay } 296*4e2b4712SSatish Balay /* backward solve the upper triangular */ 297*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 298*4e2b4712SSatish Balay v = aa + 16*diag[i] + 16; 299*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 300*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 301*4e2b4712SSatish Balay idt = 4*i; 302*4e2b4712SSatish Balay sum1 = tmp[idt]; sum2 = tmp[1+idt]; 303*4e2b4712SSatish Balay sum3 = tmp[2+idt];sum4 = tmp[3+idt]; 304*4e2b4712SSatish Balay while (nz--) { 305*4e2b4712SSatish Balay idx = 4*(*vi++); 306*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; 307*4e2b4712SSatish Balay x3 = tmp[2+idx]; x4 = tmp[3+idx]; 308*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 309*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 310*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 311*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 312*4e2b4712SSatish Balay v += 16; 313*4e2b4712SSatish Balay } 314*4e2b4712SSatish Balay idc = 4*(*c--); 315*4e2b4712SSatish Balay v = aa + 16*diag[i]; 316*4e2b4712SSatish Balay x[idc] = tmp[idt] = v[0]*sum1+v[4]*sum2+v[8]*sum3+v[12]*sum4; 317*4e2b4712SSatish Balay x[1+idc] = tmp[1+idt] = v[1]*sum1+v[5]*sum2+v[9]*sum3+v[13]*sum4; 318*4e2b4712SSatish Balay x[2+idc] = tmp[2+idt] = v[2]*sum1+v[6]*sum2+v[10]*sum3+v[14]*sum4; 319*4e2b4712SSatish Balay x[3+idc] = tmp[3+idt] = v[3]*sum1+v[7]*sum2+v[11]*sum3+v[15]*sum4; 320*4e2b4712SSatish Balay } 321*4e2b4712SSatish Balay 322*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 323*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 324*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 325*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 326*4e2b4712SSatish Balay PLogFlops(2*16*(a->nz) - a->n); 327*4e2b4712SSatish Balay PetscFunctionReturn(0); 328*4e2b4712SSatish Balay } 329*4e2b4712SSatish Balay 330*4e2b4712SSatish Balay /* 331*4e2b4712SSatish Balay Special case where the matrix was ILU(0) factored in the natural 332*4e2b4712SSatish Balay ordering. This eliminates the need for the column and row permutation. 333*4e2b4712SSatish Balay */ 334*4e2b4712SSatish Balay #undef __FUNC__ 335*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_4_NaturalOrdering" 336*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat A,Vec bb,Vec xx) 337*4e2b4712SSatish Balay { 338*4e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 339*4e2b4712SSatish Balay int i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt; 340*4e2b4712SSatish Balay int *diag = a->diag,jdx; 341*4e2b4712SSatish Balay Scalar *aa=a->a,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 342*4e2b4712SSatish Balay register Scalar *x,*b,*v; 343*4e2b4712SSatish Balay 344*4e2b4712SSatish Balay PetscFunctionBegin; 345*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 346*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 347*4e2b4712SSatish Balay 348*4e2b4712SSatish Balay /* forward solve the lower triangular */ 349*4e2b4712SSatish Balay idx = 0; 350*4e2b4712SSatish Balay x[0] = b[idx]; x[1] = b[1+idx]; x[2] = b[2+idx]; x[3] = b[3+idx]; 351*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 352*4e2b4712SSatish Balay v = aa + 16*ai[i]; 353*4e2b4712SSatish Balay vi = aj + ai[i]; 354*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 355*4e2b4712SSatish Balay idx = 4*i; 356*4e2b4712SSatish Balay sum1 = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx]; 357*4e2b4712SSatish Balay while (nz--) { 358*4e2b4712SSatish Balay jdx = 4*(*vi++); 359*4e2b4712SSatish Balay x1 = x[jdx];x2 = x[1+jdx];x3 = x[2+jdx];x4 = x[3+jdx]; 360*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 361*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 362*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 363*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 364*4e2b4712SSatish Balay v += 16; 365*4e2b4712SSatish Balay } 366*4e2b4712SSatish Balay x[idx] = sum1; 367*4e2b4712SSatish Balay x[1+idx] = sum2; 368*4e2b4712SSatish Balay x[2+idx] = sum3; 369*4e2b4712SSatish Balay x[3+idx] = sum4; 370*4e2b4712SSatish Balay } 371*4e2b4712SSatish Balay /* backward solve the upper triangular */ 372*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 373*4e2b4712SSatish Balay v = aa + 16*diag[i] + 16; 374*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 375*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 376*4e2b4712SSatish Balay idt = 4*i; 377*4e2b4712SSatish Balay sum1 = x[idt]; sum2 = x[1+idt]; 378*4e2b4712SSatish Balay sum3 = x[2+idt];sum4 = x[3+idt]; 379*4e2b4712SSatish Balay while (nz--) { 380*4e2b4712SSatish Balay idx = 4*(*vi++); 381*4e2b4712SSatish Balay x1 = x[idx]; x2 = x[1+idx];x3 = x[2+idx]; x4 = x[3+idx]; 382*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 383*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 384*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 385*4e2b4712SSatish Balay sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 386*4e2b4712SSatish Balay v += 16; 387*4e2b4712SSatish Balay } 388*4e2b4712SSatish Balay v = aa + 16*diag[i]; 389*4e2b4712SSatish Balay x[idt] = v[0]*sum1 + v[4]*sum2 + v[8]*sum3 + v[12]*sum4; 390*4e2b4712SSatish Balay x[1+idt] = v[1]*sum1 + v[5]*sum2 + v[9]*sum3 + v[13]*sum4; 391*4e2b4712SSatish Balay x[2+idt] = v[2]*sum1 + v[6]*sum2 + v[10]*sum3 + v[14]*sum4; 392*4e2b4712SSatish Balay x[3+idt] = v[3]*sum1 + v[7]*sum2 + v[11]*sum3 + v[15]*sum4; 393*4e2b4712SSatish Balay } 394*4e2b4712SSatish Balay 395*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 396*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 397*4e2b4712SSatish Balay PLogFlops(2*16*(a->nz) - a->n); 398*4e2b4712SSatish Balay PetscFunctionReturn(0); 399*4e2b4712SSatish Balay } 400*4e2b4712SSatish Balay 401*4e2b4712SSatish Balay 402*4e2b4712SSatish Balay #undef __FUNC__ 403*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_3" 404*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_3(Mat A,Vec bb,Vec xx) 405*4e2b4712SSatish Balay { 406*4e2b4712SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 407*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 408*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout; 409*4e2b4712SSatish Balay int *diag = a->diag; 410*4e2b4712SSatish Balay Scalar *aa=a->a,sum1,sum2,sum3,x1,x2,x3; 411*4e2b4712SSatish Balay register Scalar *x,*b,*tmp,*v; 412*4e2b4712SSatish Balay 413*4e2b4712SSatish Balay PetscFunctionBegin; 414*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 415*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 416*4e2b4712SSatish Balay tmp = a->solve_work; 417*4e2b4712SSatish Balay 418*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 419*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 420*4e2b4712SSatish Balay 421*4e2b4712SSatish Balay /* forward solve the lower triangular */ 422*4e2b4712SSatish Balay idx = 3*(*r++); 423*4e2b4712SSatish Balay tmp[0] = b[idx]; tmp[1] = b[1+idx]; tmp[2] = b[2+idx]; 424*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 425*4e2b4712SSatish Balay v = aa + 9*ai[i]; 426*4e2b4712SSatish Balay vi = aj + ai[i]; 427*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 428*4e2b4712SSatish Balay idx = 3*(*r++); 429*4e2b4712SSatish Balay sum1 = b[idx]; sum2 = b[1+idx]; sum3 = b[2+idx]; 430*4e2b4712SSatish Balay while (nz--) { 431*4e2b4712SSatish Balay idx = 3*(*vi++); 432*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; x3 = tmp[2+idx]; 433*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 434*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 435*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 436*4e2b4712SSatish Balay v += 9; 437*4e2b4712SSatish Balay } 438*4e2b4712SSatish Balay idx = 3*i; 439*4e2b4712SSatish Balay tmp[idx] = sum1; tmp[1+idx] = sum2; tmp[2+idx] = sum3; 440*4e2b4712SSatish Balay } 441*4e2b4712SSatish Balay /* backward solve the upper triangular */ 442*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 443*4e2b4712SSatish Balay v = aa + 9*diag[i] + 9; 444*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 445*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 446*4e2b4712SSatish Balay idt = 3*i; 447*4e2b4712SSatish Balay sum1 = tmp[idt]; sum2 = tmp[1+idt]; sum3 = tmp[2+idt]; 448*4e2b4712SSatish Balay while (nz--) { 449*4e2b4712SSatish Balay idx = 3*(*vi++); 450*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; x3 = tmp[2+idx]; 451*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 452*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 453*4e2b4712SSatish Balay sum3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 454*4e2b4712SSatish Balay v += 9; 455*4e2b4712SSatish Balay } 456*4e2b4712SSatish Balay idc = 3*(*c--); 457*4e2b4712SSatish Balay v = aa + 9*diag[i]; 458*4e2b4712SSatish Balay x[idc] = tmp[idt] = v[0]*sum1 + v[3]*sum2 + v[6]*sum3; 459*4e2b4712SSatish Balay x[1+idc] = tmp[1+idt] = v[1]*sum1 + v[4]*sum2 + v[7]*sum3; 460*4e2b4712SSatish Balay x[2+idc] = tmp[2+idt] = v[2]*sum1 + v[5]*sum2 + v[8]*sum3; 461*4e2b4712SSatish Balay } 462*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 463*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 464*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 465*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 466*4e2b4712SSatish Balay PLogFlops(2*9*(a->nz) - a->n); 467*4e2b4712SSatish Balay PetscFunctionReturn(0); 468*4e2b4712SSatish Balay } 469*4e2b4712SSatish Balay 470*4e2b4712SSatish Balay #undef __FUNC__ 471*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_2" 472*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_2(Mat A,Vec bb,Vec xx) 473*4e2b4712SSatish Balay { 474*4e2b4712SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 475*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 476*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout; 477*4e2b4712SSatish Balay int *diag = a->diag; 478*4e2b4712SSatish Balay Scalar *aa=a->a,sum1,sum2,x1,x2; 479*4e2b4712SSatish Balay register Scalar *x,*b,*tmp,*v; 480*4e2b4712SSatish Balay 481*4e2b4712SSatish Balay PetscFunctionBegin; 482*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 483*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 484*4e2b4712SSatish Balay tmp = a->solve_work; 485*4e2b4712SSatish Balay 486*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 487*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 488*4e2b4712SSatish Balay 489*4e2b4712SSatish Balay /* forward solve the lower triangular */ 490*4e2b4712SSatish Balay idx = 2*(*r++); 491*4e2b4712SSatish Balay tmp[0] = b[idx]; tmp[1] = b[1+idx]; 492*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 493*4e2b4712SSatish Balay v = aa + 4*ai[i]; 494*4e2b4712SSatish Balay vi = aj + ai[i]; 495*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 496*4e2b4712SSatish Balay idx = 2*(*r++); 497*4e2b4712SSatish Balay sum1 = b[idx]; sum2 = b[1+idx]; 498*4e2b4712SSatish Balay while (nz--) { 499*4e2b4712SSatish Balay idx = 2*(*vi++); 500*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; 501*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[2]*x2; 502*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[3]*x2; 503*4e2b4712SSatish Balay v += 4; 504*4e2b4712SSatish Balay } 505*4e2b4712SSatish Balay idx = 2*i; 506*4e2b4712SSatish Balay tmp[idx] = sum1; tmp[1+idx] = sum2; 507*4e2b4712SSatish Balay } 508*4e2b4712SSatish Balay /* backward solve the upper triangular */ 509*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 510*4e2b4712SSatish Balay v = aa + 4*diag[i] + 4; 511*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 512*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 513*4e2b4712SSatish Balay idt = 2*i; 514*4e2b4712SSatish Balay sum1 = tmp[idt]; sum2 = tmp[1+idt]; 515*4e2b4712SSatish Balay while (nz--) { 516*4e2b4712SSatish Balay idx = 2*(*vi++); 517*4e2b4712SSatish Balay x1 = tmp[idx]; x2 = tmp[1+idx]; 518*4e2b4712SSatish Balay sum1 -= v[0]*x1 + v[2]*x2; 519*4e2b4712SSatish Balay sum2 -= v[1]*x1 + v[3]*x2; 520*4e2b4712SSatish Balay v += 4; 521*4e2b4712SSatish Balay } 522*4e2b4712SSatish Balay idc = 2*(*c--); 523*4e2b4712SSatish Balay v = aa + 4*diag[i]; 524*4e2b4712SSatish Balay x[idc] = tmp[idt] = v[0]*sum1 + v[2]*sum2; 525*4e2b4712SSatish Balay x[1+idc] = tmp[1+idt] = v[1]*sum1 + v[3]*sum2; 526*4e2b4712SSatish Balay } 527*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 528*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 529*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 530*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 531*4e2b4712SSatish Balay PLogFlops(2*4*(a->nz) - a->n); 532*4e2b4712SSatish Balay PetscFunctionReturn(0); 533*4e2b4712SSatish Balay } 534*4e2b4712SSatish Balay 535*4e2b4712SSatish Balay 536*4e2b4712SSatish Balay #undef __FUNC__ 537*4e2b4712SSatish Balay #define __FUNC__ "MatSolve_SeqBAIJ_1" 538*4e2b4712SSatish Balay int MatSolve_SeqBAIJ_1(Mat A,Vec bb,Vec xx) 539*4e2b4712SSatish Balay { 540*4e2b4712SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 541*4e2b4712SSatish Balay IS iscol=a->col,isrow=a->row; 542*4e2b4712SSatish Balay int *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,*rout,*cout; 543*4e2b4712SSatish Balay int *diag = a->diag; 544*4e2b4712SSatish Balay Scalar *aa=a->a,sum1; 545*4e2b4712SSatish Balay register Scalar *x,*b,*tmp,*v; 546*4e2b4712SSatish Balay 547*4e2b4712SSatish Balay PetscFunctionBegin; 548*4e2b4712SSatish Balay if (!n) PetscFunctionReturn(0); 549*4e2b4712SSatish Balay 550*4e2b4712SSatish Balay VecGetArray_Fast(bb,b); 551*4e2b4712SSatish Balay VecGetArray_Fast(xx,x); 552*4e2b4712SSatish Balay tmp = a->solve_work; 553*4e2b4712SSatish Balay 554*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 555*4e2b4712SSatish Balay ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 556*4e2b4712SSatish Balay 557*4e2b4712SSatish Balay /* forward solve the lower triangular */ 558*4e2b4712SSatish Balay tmp[0] = b[*r++]; 559*4e2b4712SSatish Balay for ( i=1; i<n; i++ ) { 560*4e2b4712SSatish Balay v = aa + ai[i]; 561*4e2b4712SSatish Balay vi = aj + ai[i]; 562*4e2b4712SSatish Balay nz = diag[i] - ai[i]; 563*4e2b4712SSatish Balay sum1 = b[*r++]; 564*4e2b4712SSatish Balay while (nz--) { 565*4e2b4712SSatish Balay sum1 -= (*v++)*tmp[*vi++]; 566*4e2b4712SSatish Balay } 567*4e2b4712SSatish Balay tmp[i] = sum1; 568*4e2b4712SSatish Balay } 569*4e2b4712SSatish Balay /* backward solve the upper triangular */ 570*4e2b4712SSatish Balay for ( i=n-1; i>=0; i-- ){ 571*4e2b4712SSatish Balay v = aa + diag[i] + 1; 572*4e2b4712SSatish Balay vi = aj + diag[i] + 1; 573*4e2b4712SSatish Balay nz = ai[i+1] - diag[i] - 1; 574*4e2b4712SSatish Balay sum1 = tmp[i]; 575*4e2b4712SSatish Balay while (nz--) { 576*4e2b4712SSatish Balay sum1 -= (*v++)*tmp[*vi++]; 577*4e2b4712SSatish Balay } 578*4e2b4712SSatish Balay x[*c--] = tmp[i] = aa[diag[i]]*sum1; 579*4e2b4712SSatish Balay } 580*4e2b4712SSatish Balay 581*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr); 582*4e2b4712SSatish Balay ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr); 583*4e2b4712SSatish Balay VecRestoreArray_Fast(bb,b); 584*4e2b4712SSatish Balay VecRestoreArray_Fast(xx,x); 585*4e2b4712SSatish Balay PLogFlops(2*1*(a->nz) - a->n); 586*4e2b4712SSatish Balay PetscFunctionReturn(0); 587*4e2b4712SSatish Balay } 588*4e2b4712SSatish Balay 589*4e2b4712SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,Mat*); 590*4e2b4712SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 591*4e2b4712SSatish Balay /* ----------------------------------------------------------------*/ 592*4e2b4712SSatish Balay /* 593*4e2b4712SSatish Balay This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 594*4e2b4712SSatish Balay except that the data structure of Mat_SeqAIJ is slightly different. 595*4e2b4712SSatish Balay Not a good example of code reuse. 596*4e2b4712SSatish Balay */ 597*4e2b4712SSatish Balay #undef __FUNC__ 598*4e2b4712SSatish Balay #define __FUNC__ "MatILUFactorSymbolic_SeqBAIJ" 599*4e2b4712SSatish Balay int MatILUFactorSymbolic_SeqBAIJ(Mat A,IS isrow,IS iscol,double f,int levels, 600*4e2b4712SSatish Balay Mat *fact) 601*4e2b4712SSatish Balay { 602*4e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data, *b; 603*4e2b4712SSatish Balay IS isicol; 604*4e2b4712SSatish Balay int *r,*ic, ierr, prow, n = a->mbs, *ai = a->i, *aj = a->j; 605*4e2b4712SSatish Balay int *ainew,*ajnew, jmax,*fill, *xi, nz, *im,*ajfill,*flev; 606*4e2b4712SSatish Balay int *dloc, idx, row,m,fm, nzf, nzi,len, realloc = 0; 607*4e2b4712SSatish Balay int incrlev,nnz,i,bs = a->bs,bs2 = a->bs2; 608*4e2b4712SSatish Balay PetscTruth col_identity, row_identity; 609*4e2b4712SSatish Balay 610*4e2b4712SSatish Balay PetscFunctionBegin; 611*4e2b4712SSatish Balay /* special case that simply copies fill pattern */ 612*4e2b4712SSatish Balay PetscValidHeaderSpecific(isrow,IS_COOKIE); 613*4e2b4712SSatish Balay PetscValidHeaderSpecific(iscol,IS_COOKIE); 614*4e2b4712SSatish Balay ISIdentity(isrow,&row_identity); ISIdentity(iscol,&col_identity); 615*4e2b4712SSatish Balay if (levels == 0 && row_identity && col_identity) { 616*4e2b4712SSatish Balay ierr = MatConvertSameType_SeqBAIJ(A,fact,DO_NOT_COPY_VALUES); CHKERRQ(ierr); 617*4e2b4712SSatish Balay (*fact)->factor = FACTOR_LU; 618*4e2b4712SSatish Balay b = (Mat_SeqBAIJ *) (*fact)->data; 619*4e2b4712SSatish Balay if (!b->diag) { 620*4e2b4712SSatish Balay ierr = MatMarkDiag_SeqBAIJ(*fact); CHKERRQ(ierr); 621*4e2b4712SSatish Balay } 622*4e2b4712SSatish Balay b->row = isrow; 623*4e2b4712SSatish Balay b->col = iscol; 624*4e2b4712SSatish Balay b->solve_work = (Scalar *) PetscMalloc((b->m+1+b->bs)*sizeof(Scalar));CHKPTRQ(b->solve_work); 625*4e2b4712SSatish Balay /* 626*4e2b4712SSatish Balay Blocksize 4 has a special faster solver for ILU(0) factorization 627*4e2b4712SSatish Balay with natural ordering 628*4e2b4712SSatish Balay */ 629*4e2b4712SSatish Balay if (b->bs == 4) { 630*4e2b4712SSatish Balay (*fact)->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 631*4e2b4712SSatish Balay (*fact)->ops.solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 632*4e2b4712SSatish Balay } 633*4e2b4712SSatish Balay PetscFunctionReturn(0); 634*4e2b4712SSatish Balay } 635*4e2b4712SSatish Balay 636*4e2b4712SSatish Balay ierr = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr); 637*4e2b4712SSatish Balay ierr = ISGetIndices(isrow,&r); CHKERRQ(ierr); 638*4e2b4712SSatish Balay ierr = ISGetIndices(isicol,&ic); CHKERRQ(ierr); 639*4e2b4712SSatish Balay 640*4e2b4712SSatish Balay /* get new row pointers */ 641*4e2b4712SSatish Balay ainew = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(ainew); 642*4e2b4712SSatish Balay ainew[0] = 0; 643*4e2b4712SSatish Balay /* don't know how many column pointers are needed so estimate */ 644*4e2b4712SSatish Balay jmax = (int) (f*ai[n] + 1); 645*4e2b4712SSatish Balay ajnew = (int *) PetscMalloc( (jmax)*sizeof(int) ); CHKPTRQ(ajnew); 646*4e2b4712SSatish Balay /* ajfill is level of fill for each fill entry */ 647*4e2b4712SSatish Balay ajfill = (int *) PetscMalloc( (jmax)*sizeof(int) ); CHKPTRQ(ajfill); 648*4e2b4712SSatish Balay /* fill is a linked list of nonzeros in active row */ 649*4e2b4712SSatish Balay fill = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(fill); 650*4e2b4712SSatish Balay /* im is level for each filled value */ 651*4e2b4712SSatish Balay im = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(im); 652*4e2b4712SSatish Balay /* dloc is location of diagonal in factor */ 653*4e2b4712SSatish Balay dloc = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(dloc); 654*4e2b4712SSatish Balay dloc[0] = 0; 655*4e2b4712SSatish Balay for ( prow=0; prow<n; prow++ ) { 656*4e2b4712SSatish Balay /* first copy previous fill into linked list */ 657*4e2b4712SSatish Balay nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 658*4e2b4712SSatish Balay if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,1,"Empty row in matrix"); 659*4e2b4712SSatish Balay xi = aj + ai[r[prow]]; 660*4e2b4712SSatish Balay fill[n] = n; 661*4e2b4712SSatish Balay while (nz--) { 662*4e2b4712SSatish Balay fm = n; 663*4e2b4712SSatish Balay idx = ic[*xi++]; 664*4e2b4712SSatish Balay do { 665*4e2b4712SSatish Balay m = fm; 666*4e2b4712SSatish Balay fm = fill[m]; 667*4e2b4712SSatish Balay } while (fm < idx); 668*4e2b4712SSatish Balay fill[m] = idx; 669*4e2b4712SSatish Balay fill[idx] = fm; 670*4e2b4712SSatish Balay im[idx] = 0; 671*4e2b4712SSatish Balay } 672*4e2b4712SSatish Balay nzi = 0; 673*4e2b4712SSatish Balay row = fill[n]; 674*4e2b4712SSatish Balay while ( row < prow ) { 675*4e2b4712SSatish Balay incrlev = im[row] + 1; 676*4e2b4712SSatish Balay nz = dloc[row]; 677*4e2b4712SSatish Balay xi = ajnew + ainew[row] + nz; 678*4e2b4712SSatish Balay flev = ajfill + ainew[row] + nz + 1; 679*4e2b4712SSatish Balay nnz = ainew[row+1] - ainew[row] - nz - 1; 680*4e2b4712SSatish Balay if (*xi++ != row) { 681*4e2b4712SSatish Balay SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Zero pivot: try running with -pc_ilu_nonzeros_along_diagonal"); 682*4e2b4712SSatish Balay } 683*4e2b4712SSatish Balay fm = row; 684*4e2b4712SSatish Balay while (nnz-- > 0) { 685*4e2b4712SSatish Balay idx = *xi++; 686*4e2b4712SSatish Balay if (*flev + incrlev > levels) { 687*4e2b4712SSatish Balay flev++; 688*4e2b4712SSatish Balay continue; 689*4e2b4712SSatish Balay } 690*4e2b4712SSatish Balay do { 691*4e2b4712SSatish Balay m = fm; 692*4e2b4712SSatish Balay fm = fill[m]; 693*4e2b4712SSatish Balay } while (fm < idx); 694*4e2b4712SSatish Balay if (fm != idx) { 695*4e2b4712SSatish Balay im[idx] = *flev + incrlev; 696*4e2b4712SSatish Balay fill[m] = idx; 697*4e2b4712SSatish Balay fill[idx] = fm; 698*4e2b4712SSatish Balay fm = idx; 699*4e2b4712SSatish Balay nzf++; 700*4e2b4712SSatish Balay } 701*4e2b4712SSatish Balay else { 702*4e2b4712SSatish Balay if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 703*4e2b4712SSatish Balay } 704*4e2b4712SSatish Balay flev++; 705*4e2b4712SSatish Balay } 706*4e2b4712SSatish Balay row = fill[row]; 707*4e2b4712SSatish Balay nzi++; 708*4e2b4712SSatish Balay } 709*4e2b4712SSatish Balay /* copy new filled row into permanent storage */ 710*4e2b4712SSatish Balay ainew[prow+1] = ainew[prow] + nzf; 711*4e2b4712SSatish Balay if (ainew[prow+1] > jmax) { 712*4e2b4712SSatish Balay /* allocate a longer ajnew */ 713*4e2b4712SSatish Balay int maxadd; 714*4e2b4712SSatish Balay maxadd = (int) (((f*ai[n]+1)*(n-prow+5))/n); 715*4e2b4712SSatish Balay if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 716*4e2b4712SSatish Balay jmax += maxadd; 717*4e2b4712SSatish Balay xi = (int *) PetscMalloc( jmax*sizeof(int) );CHKPTRQ(xi); 718*4e2b4712SSatish Balay PetscMemcpy(xi,ajnew,ainew[prow]*sizeof(int)); 719*4e2b4712SSatish Balay PetscFree(ajnew); 720*4e2b4712SSatish Balay ajnew = xi; 721*4e2b4712SSatish Balay /* allocate a longer ajfill */ 722*4e2b4712SSatish Balay xi = (int *) PetscMalloc( jmax*sizeof(int) );CHKPTRQ(xi); 723*4e2b4712SSatish Balay PetscMemcpy(xi,ajfill,ainew[prow]*sizeof(int)); 724*4e2b4712SSatish Balay PetscFree(ajfill); 725*4e2b4712SSatish Balay ajfill = xi; 726*4e2b4712SSatish Balay realloc++; 727*4e2b4712SSatish Balay } 728*4e2b4712SSatish Balay xi = ajnew + ainew[prow]; 729*4e2b4712SSatish Balay flev = ajfill + ainew[prow]; 730*4e2b4712SSatish Balay dloc[prow] = nzi; 731*4e2b4712SSatish Balay fm = fill[n]; 732*4e2b4712SSatish Balay while (nzf--) { 733*4e2b4712SSatish Balay *xi++ = fm; 734*4e2b4712SSatish Balay *flev++ = im[fm]; 735*4e2b4712SSatish Balay fm = fill[fm]; 736*4e2b4712SSatish Balay } 737*4e2b4712SSatish Balay } 738*4e2b4712SSatish Balay PetscFree(ajfill); 739*4e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr); 740*4e2b4712SSatish Balay ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr); 741*4e2b4712SSatish Balay ierr = ISDestroy(isicol); CHKERRQ(ierr); 742*4e2b4712SSatish Balay PetscFree(fill); PetscFree(im); 743*4e2b4712SSatish Balay 744*4e2b4712SSatish Balay { 745*4e2b4712SSatish Balay double af = ((double)ainew[n])/((double)ai[n]); 746*4e2b4712SSatish Balay PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:Reallocs %d Fill ratio:given %g needed %g\n", 747*4e2b4712SSatish Balay realloc,f,af); 748*4e2b4712SSatish Balay PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:Run with -pc_ilu_fill %g or use \n",af); 749*4e2b4712SSatish Balay PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:PCILUSetFill(pc,%g);\n",af); 750*4e2b4712SSatish Balay PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:for best performance.\n"); 751*4e2b4712SSatish Balay } 752*4e2b4712SSatish Balay 753*4e2b4712SSatish Balay /* put together the new matrix */ 754*4e2b4712SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,bs*n,bs*n,0,PETSC_NULL,fact);CHKERRQ(ierr); 755*4e2b4712SSatish Balay b = (Mat_SeqBAIJ *) (*fact)->data; 756*4e2b4712SSatish Balay PetscFree(b->imax); 757*4e2b4712SSatish Balay b->singlemalloc = 0; 758*4e2b4712SSatish Balay len = bs2*ainew[n]*sizeof(Scalar); 759*4e2b4712SSatish Balay /* the next line frees the default space generated by the Create() */ 760*4e2b4712SSatish Balay PetscFree(b->a); PetscFree(b->ilen); 761*4e2b4712SSatish Balay b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 762*4e2b4712SSatish Balay b->j = ajnew; 763*4e2b4712SSatish Balay b->i = ainew; 764*4e2b4712SSatish Balay for ( i=0; i<n; i++ ) dloc[i] += ainew[i]; 765*4e2b4712SSatish Balay b->diag = dloc; 766*4e2b4712SSatish Balay b->ilen = 0; 767*4e2b4712SSatish Balay b->imax = 0; 768*4e2b4712SSatish Balay b->row = isrow; 769*4e2b4712SSatish Balay b->col = iscol; 770*4e2b4712SSatish Balay b->solve_work = (Scalar *) PetscMalloc( (bs*n+bs)*sizeof(Scalar));CHKPTRQ(b->solve_work); 771*4e2b4712SSatish Balay /* In b structure: Free imax, ilen, old a, old j. 772*4e2b4712SSatish Balay Allocate dloc, solve_work, new a, new j */ 773*4e2b4712SSatish Balay PLogObjectMemory(*fact,(ainew[n]-n)*(sizeof(int))+bs2*ainew[n]*sizeof(Scalar)); 774*4e2b4712SSatish Balay b->maxnz = b->nz = ainew[n]; 775*4e2b4712SSatish Balay (*fact)->factor = FACTOR_LU; 776*4e2b4712SSatish Balay 777*4e2b4712SSatish Balay (*fact)->info.factor_mallocs = realloc; 778*4e2b4712SSatish Balay (*fact)->info.fill_ratio_given = f; 779*4e2b4712SSatish Balay (*fact)->info.fill_ratio_needed = ((double)ainew[n])/((double)ai[prow]); 780*4e2b4712SSatish Balay 781*4e2b4712SSatish Balay PetscFunctionReturn(0); 782*4e2b4712SSatish Balay } 783*4e2b4712SSatish Balay 784*4e2b4712SSatish Balay 785*4e2b4712SSatish Balay 786*4e2b4712SSatish Balay 787