PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; PetscErrorCode ierr; PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; PetscInt *ajtmpold,*ajtmp,nz,row; PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15; MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25; MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15; MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25; MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15; MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36; MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36; MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36; MatScalar *ba = b->a,*aa = a->a; PetscReal shift = info->shiftamount; PetscFunctionBegin; ierr = PetscMalloc(36*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); for (i=0; i<n; i++) { nz = bi[i+1] - bi[i]; ajtmp = bj + bi[i]; for (j=0; j<nz; j++) { x = rtmp+36*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0; x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0 ; x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0 ; x[34] = x[35] = 0.0 ; } /* load in initial (unfactored row) */ nz = ai[i+1] - ai[i]; ajtmpold = aj + ai[i]; v = aa + 36*ai[i]; for (j=0; j<nz; j++) { x = rtmp+36*ajtmpold[j]; x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23]; x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27]; x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31]; x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35]; v += 36; } row = *ajtmp++; while (row < i) { pc = rtmp + 36*row; p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23]; p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27]; p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31]; p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35]; if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 || p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 || p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) { pv = ba + 36*diag_offset[row]; pj = bj + diag_offset[row] + 1; x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27]; x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31]; x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35]; pc[0] = m1 = p1*x1 + p7*x2 + p13*x3 + p19*x4 + p25*x5 + p31*x6; pc[1] = m2 = p2*x1 + p8*x2 + p14*x3 + p20*x4 + p26*x5 + p32*x6; pc[2] = m3 = p3*x1 + p9*x2 + p15*x3 + p21*x4 + p27*x5 + p33*x6; pc[3] = m4 = p4*x1 + p10*x2 + p16*x3 + p22*x4 + p28*x5 + p34*x6; pc[4] = m5 = p5*x1 + p11*x2 + p17*x3 + p23*x4 + p29*x5 + p35*x6; pc[5] = m6 = p6*x1 + p12*x2 + p18*x3 + p24*x4 + p30*x5 + p36*x6; pc[6] = m7 = p1*x7 + p7*x8 + p13*x9 + p19*x10 + p25*x11 + p31*x12; pc[7] = m8 = p2*x7 + p8*x8 + p14*x9 + p20*x10 + p26*x11 + p32*x12; pc[8] = m9 = p3*x7 + p9*x8 + p15*x9 + p21*x10 + p27*x11 + p33*x12; pc[9] = m10 = p4*x7 + p10*x8 + p16*x9 + p22*x10 + p28*x11 + p34*x12; pc[10] = m11 = p5*x7 + p11*x8 + p17*x9 + p23*x10 + p29*x11 + p35*x12; pc[11] = m12 = p6*x7 + p12*x8 + p18*x9 + p24*x10 + p30*x11 + p36*x12; pc[12] = m13 = p1*x13 + p7*x14 + p13*x15 + p19*x16 + p25*x17 + p31*x18; pc[13] = m14 = p2*x13 + p8*x14 + p14*x15 + p20*x16 + p26*x17 + p32*x18; pc[14] = m15 = p3*x13 + p9*x14 + p15*x15 + p21*x16 + p27*x17 + p33*x18; pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18; pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18; pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18; pc[18] = m19 = p1*x19 + p7*x20 + p13*x21 + p19*x22 + p25*x23 + p31*x24; pc[19] = m20 = p2*x19 + p8*x20 + p14*x21 + p20*x22 + p26*x23 + p32*x24; pc[20] = m21 = p3*x19 + p9*x20 + p15*x21 + p21*x22 + p27*x23 + p33*x24; pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24; pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24; pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24; pc[24] = m25 = p1*x25 + p7*x26 + p13*x27 + p19*x28 + p25*x29 + p31*x30; pc[25] = m26 = p2*x25 + p8*x26 + p14*x27 + p20*x28 + p26*x29 + p32*x30; pc[26] = m27 = p3*x25 + p9*x26 + p15*x27 + p21*x28 + p27*x29 + p33*x30; pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30; pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30; pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30; pc[30] = m31 = p1*x31 + p7*x32 + p13*x33 + p19*x34 + p25*x35 + p31*x36; pc[31] = m32 = p2*x31 + p8*x32 + p14*x33 + p20*x34 + p26*x35 + p32*x36; pc[32] = m33 = p3*x31 + p9*x32 + p15*x33 + p21*x34 + p27*x35 + p33*x36; pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36; pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36; pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36; nz = bi[row+1] - diag_offset[row] - 1; pv += 36; for (j=0; j<nz; j++) { x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27]; x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31]; x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35]; x = rtmp + 36*pj[j]; x[0] -= m1*x1 + m7*x2 + m13*x3 + m19*x4 + m25*x5 + m31*x6; x[1] -= m2*x1 + m8*x2 + m14*x3 + m20*x4 + m26*x5 + m32*x6; x[2] -= m3*x1 + m9*x2 + m15*x3 + m21*x4 + m27*x5 + m33*x6; x[3] -= m4*x1 + m10*x2 + m16*x3 + m22*x4 + m28*x5 + m34*x6; x[4] -= m5*x1 + m11*x2 + m17*x3 + m23*x4 + m29*x5 + m35*x6; x[5] -= m6*x1 + m12*x2 + m18*x3 + m24*x4 + m30*x5 + m36*x6; x[6] -= m1*x7 + m7*x8 + m13*x9 + m19*x10 + m25*x11 + m31*x12; x[7] -= m2*x7 + m8*x8 + m14*x9 + m20*x10 + m26*x11 + m32*x12; x[8] -= m3*x7 + m9*x8 + m15*x9 + m21*x10 + m27*x11 + m33*x12; x[9] -= m4*x7 + m10*x8 + m16*x9 + m22*x10 + m28*x11 + m34*x12; x[10] -= m5*x7 + m11*x8 + m17*x9 + m23*x10 + m29*x11 + m35*x12; x[11] -= m6*x7 + m12*x8 + m18*x9 + m24*x10 + m30*x11 + m36*x12; x[12] -= m1*x13 + m7*x14 + m13*x15 + m19*x16 + m25*x17 + m31*x18; x[13] -= m2*x13 + m8*x14 + m14*x15 + m20*x16 + m26*x17 + m32*x18; x[14] -= m3*x13 + m9*x14 + m15*x15 + m21*x16 + m27*x17 + m33*x18; x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18; x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18; x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18; x[18] -= m1*x19 + m7*x20 + m13*x21 + m19*x22 + m25*x23 + m31*x24; x[19] -= m2*x19 + m8*x20 + m14*x21 + m20*x22 + m26*x23 + m32*x24; x[20] -= m3*x19 + m9*x20 + m15*x21 + m21*x22 + m27*x23 + m33*x24; x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24; x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24; x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24; x[24] -= m1*x25 + m7*x26 + m13*x27 + m19*x28 + m25*x29 + m31*x30; x[25] -= m2*x25 + m8*x26 + m14*x27 + m20*x28 + m26*x29 + m32*x30; x[26] -= m3*x25 + m9*x26 + m15*x27 + m21*x28 + m27*x29 + m33*x30; x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30; x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30; x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30; x[30] -= m1*x31 + m7*x32 + m13*x33 + m19*x34 + m25*x35 + m31*x36; x[31] -= m2*x31 + m8*x32 + m14*x33 + m20*x34 + m26*x35 + m32*x36; x[32] -= m3*x31 + m9*x32 + m15*x33 + m21*x34 + m27*x35 + m33*x36; x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36; x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36; x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36; pv += 36; } ierr = PetscLogFlops(432.0*nz+396.0);CHKERRQ(ierr); } row = *ajtmp++; } /* finished row so stick it into b->a */ pv = ba + 36*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; j<nz; j++) { x = rtmp+36*pj[j]; pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27]; pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31]; pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35]; pv += 36; } /* invert diagonal block */ w = ba + 36*diag_offset[i]; ierr = PetscKernel_A_gets_inverse_A_6(w,shift);CHKERRQ(ierr); } ierr = PetscFree(rtmp);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace; C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace; C->assembled = PETSC_TRUE; ierr = PetscLogFlops(1.333333333333*6*6*6*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); }
PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) { Mat C=B; Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; PetscErrorCode ierr; PetscInt i,j,k,nz,nzL,row; const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; PetscInt flg; PetscReal shift = info->shiftamount; PetscFunctionBegin; /* generate work space needed by the factorization */ ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); for (i=0; i<n; i++){ /* zero rtmp */ /* L part */ nz = bi[i+1] - bi[i]; bjtmp = bj + bi[i]; for (j=0; j<nz; j++){ ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); } /* U part */ nz = bdiag[i] - bdiag[i+1]; bjtmp = bj + bdiag[i+1]+1; for (j=0; j<nz; j++){ ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); } /* load in initial (unfactored row) */ nz = ai[i+1] - ai[i]; ajtmp = aj + ai[i]; v = aa + bs2*ai[i]; for (j=0; j<nz; j++) { ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); } /* elimination */ bjtmp = bj + bi[i]; nzL = bi[i+1] - bi[i]; for(k=0;k < nzL;k++) { row = bjtmp[k]; pc = rtmp + bs2*row; for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} if (flg) { pv = b->a + bs2*bdiag[row]; /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ ierr = PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);CHKERRQ(ierr); pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ pv = b->a + bs2*(bdiag[row+1]+1); nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ for (j=0; j<nz; j++) { /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ v = rtmp + bs2*pj[j]; ierr = PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);CHKERRQ(ierr); pv += bs2; } ierr = PetscLogFlops(432*nz+396);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ } } /* finished row so stick it into b->a */ /* L part */ pv = b->a + bs2*bi[i] ; pj = b->j + bi[i] ; nz = bi[i+1] - bi[i]; for (j=0; j<nz; j++) { ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); } /* Mark diagonal and invert diagonal for simplier triangular solves */ pv = b->a + bs2*bdiag[i]; pj = b->j + bdiag[i]; ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ ierr = PetscKernel_A_gets_inverse_A_6(pv,shift);CHKERRQ(ierr); /* U part */ pv = b->a + bs2*(bdiag[i+1]+1); pj = b->j + bdiag[i+1]+1; nz = bdiag[i] - bdiag[i+1] - 1; for (j=0; j<nz; j++){ ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); } } ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering; C->assembled = PETSC_TRUE; ierr = PetscLogFlops(1.333333333333*6*6*6*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); }
PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; PetscErrorCode ierr; PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*d,*w,*wp,u0,u1,u2,u3,u4,u5,u6,u7,u8,u9,u10,u11,u12; MatScalar u13,u14,u15,u16,u17,u18,u19,u20,u21,u22,u23,u24,u25,u26,u27; MatScalar u28,u29,u30,u31,u32,u33,u34,u35; MatScalar d0,d1,d2,d3,d4,d5,d6,d7,d8,d9,d10,d11,d12; MatScalar d13,d14,d15,d16,d17,d18,d19,d20,d21,d22,d23,d24,d25,d26,d27; MatScalar d28,d29,d30,d31,d32,d33,d34,d35; MatScalar m0,m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12; MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25,m26,m27; MatScalar m28,m29,m30,m31,m32,m33,m34,m35; PetscReal shift = info->shiftamount; PetscBool allowzeropivot,zeropivotdetected; PetscFunctionBegin; /* initialization */ allowzeropivot = PetscNot(A->erroriffailure); ierr = PetscCalloc1(36*mbs,&w);CHKERRQ(ierr); ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); for (i=0; i<mbs; i++) { jl[i] = mbs; il[0] = 0; } ierr = PetscMalloc2(36,&dk,36,&uik);CHKERRQ(ierr); ai = a->i; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; k<mbs; k++) { /*initialize k-th row with elements nonzero in row k of A */ jmin = ai[k]; jmax = ai[k+1]; if (jmin < jmax) { ap = aa + jmin*36; for (j = jmin; j < jmax; j++) { vj = aj[j]; /* block col. index */ wp = w + vj*36; for (i=0; i<36; i++) *wp++ = *ap++; } } /* modify k-th row by adding in those rows i with U(i,k) != 0 */ ierr = PetscMemcpy(dk,w+k*36,36*sizeof(MatScalar));CHKERRQ(ierr); i = jl[k]; /* first row to be added to k_th row */ while (i < mbs) { nexti = jl[i]; /* next row to be added to k_th row */ /* compute multiplier */ ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ /* uik = -inv(Di)*U_bar(i,k) */ d = ba + i*36; u = ba + ili*36; u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6]; u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13]; u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20]; u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27]; u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34]; u35 = u[35]; d0 = d[0]; d1 = d[1]; d2 = d[2]; d3 = d[3]; d4 = d[4]; d5 = d[5]; d6 = d[6]; d7 = d[7]; d8 = d[8]; d9 = d[9]; d10 = d[10]; d11 = d[11]; d12 = d[12]; d13 = d[13]; d14 = d[14]; d15 = d[15]; d16 = d[16]; d17 = d[17]; d18 = d[18]; d19 = d[19]; d20 = d[20]; d21 = d[21]; d22 = d[22]; d23 = d[23]; d24 = d[24]; d25 = d[25]; d26 = d[26]; d27 = d[27]; d28 = d[28]; d29 = d[29]; d30 = d[30]; d31 = d[31]; d32 = d[32]; d33 = d[33]; d34 = d[34]; d35 = d[35]; m0 = uik[0] = -(d0*u0 + d6*u1 + d12*u2 + d18*u3 + d24*u4 + d30*u5); m1 = uik[1] = -(d1*u0 + d7*u1 + d13*u2 + d19*u3 + d25*u4 + d31*u5); m2 = uik[2] = -(d2*u0 + d8*u1 + d14*u2 + d20*u3 + d26*u4 + d32*u5); m3 = uik[3] = -(d3*u0 + d9*u1 + d15*u2 + d21*u3 + d27*u4 + d33*u5); m4 = uik[4] = -(d4*u0+ d10*u1 + d16*u2 + d22*u3 + d28*u4 + d34*u5); m5 = uik[5] = -(d5*u0+ d11*u1 + d17*u2 + d23*u3 + d29*u4 + d35*u5); m6 = uik[6] = -(d0*u6 + d6*u7 + d12*u8 + d18*u9 + d24*u10 + d30*u11); m7 = uik[7] = -(d1*u6 + d7*u7 + d13*u8 + d19*u9 + d25*u10 + d31*u11); m8 = uik[8] = -(d2*u6 + d8*u7 + d14*u8 + d20*u9 + d26*u10 + d32*u11); m9 = uik[9] = -(d3*u6 + d9*u7 + d15*u8 + d21*u9 + d27*u10 + d33*u11); m10 = uik[10]= -(d4*u6+ d10*u7 + d16*u8 + d22*u9 + d28*u10 + d34*u11); m11 = uik[11]= -(d5*u6+ d11*u7 + d17*u8 + d23*u9 + d29*u10 + d35*u11); m12 = uik[12] = -(d0*u12 + d6*u13 + d12*u14 + d18*u15 + d24*u16 + d30*u17); m13 = uik[13] = -(d1*u12 + d7*u13 + d13*u14 + d19*u15 + d25*u16 + d31*u17); m14 = uik[14] = -(d2*u12 + d8*u13 + d14*u14 + d20*u15 + d26*u16 + d32*u17); m15 = uik[15] = -(d3*u12 + d9*u13 + d15*u14 + d21*u15 + d27*u16 + d33*u17); m16 = uik[16] = -(d4*u12+ d10*u13 + d16*u14 + d22*u15 + d28*u16 + d34*u17); m17 = uik[17] = -(d5*u12+ d11*u13 + d17*u14 + d23*u15 + d29*u16 + d35*u17); m18 = uik[18] = -(d0*u18 + d6*u19 + d12*u20 + d18*u21 + d24*u22 + d30*u23); m19 = uik[19] = -(d1*u18 + d7*u19 + d13*u20 + d19*u21 + d25*u22 + d31*u23); m20 = uik[20] = -(d2*u18 + d8*u19 + d14*u20 + d20*u21 + d26*u22 + d32*u23); m21 = uik[21] = -(d3*u18 + d9*u19 + d15*u20 + d21*u21 + d27*u22 + d33*u23); m22 = uik[22] = -(d4*u18+ d10*u19 + d16*u20 + d22*u21 + d28*u22 + d34*u23); m23 = uik[23] = -(d5*u18+ d11*u19 + d17*u20 + d23*u21 + d29*u22 + d35*u23); m24 = uik[24] = -(d0*u24 + d6*u25 + d12*u26 + d18*u27 + d24*u28 + d30*u29); m25 = uik[25] = -(d1*u24 + d7*u25 + d13*u26 + d19*u27 + d25*u28 + d31*u29); m26 = uik[26] = -(d2*u24 + d8*u25 + d14*u26 + d20*u27 + d26*u28 + d32*u29); m27 = uik[27] = -(d3*u24 + d9*u25 + d15*u26 + d21*u27 + d27*u28 + d33*u29); m28 = uik[28] = -(d4*u24+ d10*u25 + d16*u26 + d22*u27 + d28*u28 + d34*u29); m29 = uik[29] = -(d5*u24+ d11*u25 + d17*u26 + d23*u27 + d29*u28 + d35*u29); m30 = uik[30] = -(d0*u30 + d6*u31 + d12*u32 + d18*u33 + d24*u34 + d30*u35); m31 = uik[31] = -(d1*u30 + d7*u31 + d13*u32 + d19*u33 + d25*u34 + d31*u35); m32 = uik[32] = -(d2*u30 + d8*u31 + d14*u32 + d20*u33 + d26*u34 + d32*u35); m33 = uik[33] = -(d3*u30 + d9*u31 + d15*u32 + d21*u33 + d27*u34 + d33*u35); m34 = uik[34] = -(d4*u30+ d10*u31 + d16*u32 + d22*u33 + d28*u34 + d34*u35); m35 = uik[35] = -(d5*u30+ d11*u31 + d17*u32 + d23*u33 + d29*u34 + d35*u35); /* update D(k) += -U(i,k)^T * U_bar(i,k) */ dk[0] += m0*u0 + m1*u1 + m2*u2 + m3*u3 + m4*u4 + m5*u5; dk[1] += m6*u0 + m7*u1 + m8*u2 + m9*u3+ m10*u4+ m11*u5; dk[2] += m12*u0+ m13*u1+ m14*u2+ m15*u3+ m16*u4+ m17*u5; dk[3] += m18*u0+ m19*u1+ m20*u2+ m21*u3+ m22*u4+ m23*u5; dk[4] += m24*u0+ m25*u1+ m26*u2+ m27*u3+ m28*u4+ m29*u5; dk[5] += m30*u0+ m31*u1+ m32*u2+ m33*u3+ m34*u4+ m35*u5; dk[6] += m0*u6 + m1*u7 + m2*u8 + m3*u9 + m4*u10 + m5*u11; dk[7] += m6*u6 + m7*u7 + m8*u8 + m9*u9+ m10*u10+ m11*u11; dk[8] += m12*u6+ m13*u7+ m14*u8+ m15*u9+ m16*u10+ m17*u11; dk[9] += m18*u6+ m19*u7+ m20*u8+ m21*u9+ m22*u10+ m23*u11; dk[10]+= m24*u6+ m25*u7+ m26*u8+ m27*u9+ m28*u10+ m29*u11; dk[11]+= m30*u6+ m31*u7+ m32*u8+ m33*u9+ m34*u10+ m35*u11; dk[12]+= m0*u12 + m1*u13 + m2*u14 + m3*u15 + m4*u16 + m5*u17; dk[13]+= m6*u12 + m7*u13 + m8*u14 + m9*u15+ m10*u16+ m11*u17; dk[14]+= m12*u12+ m13*u13+ m14*u14+ m15*u15+ m16*u16+ m17*u17; dk[15]+= m18*u12+ m19*u13+ m20*u14+ m21*u15+ m22*u16+ m23*u17; dk[16]+= m24*u12+ m25*u13+ m26*u14+ m27*u15+ m28*u16+ m29*u17; dk[17]+= m30*u12+ m31*u13+ m32*u14+ m33*u15+ m34*u16+ m35*u17; dk[18]+= m0*u18 + m1*u19 + m2*u20 + m3*u21 + m4*u22 + m5*u23; dk[19]+= m6*u18 + m7*u19 + m8*u20 + m9*u21+ m10*u22+ m11*u23; dk[20]+= m12*u18+ m13*u19+ m14*u20+ m15*u21+ m16*u22+ m17*u23; dk[21]+= m18*u18+ m19*u19+ m20*u20+ m21*u21+ m22*u22+ m23*u23; dk[22]+= m24*u18+ m25*u19+ m26*u20+ m27*u21+ m28*u22+ m29*u23; dk[23]+= m30*u18+ m31*u19+ m32*u20+ m33*u21+ m34*u22+ m35*u23; dk[24]+= m0*u24 + m1*u25 + m2*u26 + m3*u27 + m4*u28 + m5*u29; dk[25]+= m6*u24 + m7*u25 + m8*u26 + m9*u27+ m10*u28+ m11*u29; dk[26]+= m12*u24+ m13*u25+ m14*u26+ m15*u27+ m16*u28+ m17*u29; dk[27]+= m18*u24+ m19*u25+ m20*u26+ m21*u27+ m22*u28+ m23*u29; dk[28]+= m24*u24+ m25*u25+ m26*u26+ m27*u27+ m28*u28+ m29*u29; dk[29]+= m30*u24+ m31*u25+ m32*u26+ m33*u27+ m34*u28+ m35*u29; dk[30]+= m0*u30 + m1*u31 + m2*u32 + m3*u33 + m4*u34 + m5*u35; dk[31]+= m6*u30 + m7*u31 + m8*u32 + m9*u33+ m10*u34+ m11*u35; dk[32]+= m12*u30+ m13*u31+ m14*u32+ m15*u33+ m16*u34+ m17*u35; dk[33]+= m18*u30+ m19*u31+ m20*u32+ m21*u33+ m22*u34+ m23*u35; dk[34]+= m24*u30+ m25*u31+ m26*u32+ m27*u33+ m28*u34+ m29*u35; dk[35]+= m30*u30+ m31*u31+ m32*u32+ m33*u33+ m34*u34+ m35*u35; ierr = PetscLogFlops(216.0*4.0);CHKERRQ(ierr); /* update -U(i,k) */ ierr = PetscMemcpy(ba+ili*36,uik,36*sizeof(MatScalar));CHKERRQ(ierr); /* add multiple of row i to k-th row ... */ jmin = ili + 1; jmax = bi[i+1]; if (jmin < jmax) { for (j=jmin; j<jmax; j++) { /* w += -U(i,k)^T * U_bar(i,j) */ wp = w + bj[j]*36; u = ba + j*36; u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6]; u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13]; u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20]; u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27]; u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34]; u35 = u[35]; wp[0] += m0*u0 + m1*u1 + m2*u2 + m3*u3 + m4*u4 + m5*u5; wp[1] += m6*u0 + m7*u1 + m8*u2 + m9*u3+ m10*u4+ m11*u5; wp[2] += m12*u0+ m13*u1+ m14*u2+ m15*u3+ m16*u4+ m17*u5; wp[3] += m18*u0+ m19*u1+ m20*u2+ m21*u3+ m22*u4+ m23*u5; wp[4] += m24*u0+ m25*u1+ m26*u2+ m27*u3+ m28*u4+ m29*u5; wp[5] += m30*u0+ m31*u1+ m32*u2+ m33*u3+ m34*u4+ m35*u5; wp[6] += m0*u6 + m1*u7 + m2*u8 + m3*u9 + m4*u10 + m5*u11; wp[7] += m6*u6 + m7*u7 + m8*u8 + m9*u9+ m10*u10+ m11*u11; wp[8] += m12*u6+ m13*u7+ m14*u8+ m15*u9+ m16*u10+ m17*u11; wp[9] += m18*u6+ m19*u7+ m20*u8+ m21*u9+ m22*u10+ m23*u11; wp[10]+= m24*u6+ m25*u7+ m26*u8+ m27*u9+ m28*u10+ m29*u11; wp[11]+= m30*u6+ m31*u7+ m32*u8+ m33*u9+ m34*u10+ m35*u11; wp[12]+= m0*u12 + m1*u13 + m2*u14 + m3*u15 + m4*u16 + m5*u17; wp[13]+= m6*u12 + m7*u13 + m8*u14 + m9*u15+ m10*u16+ m11*u17; wp[14]+= m12*u12+ m13*u13+ m14*u14+ m15*u15+ m16*u16+ m17*u17; wp[15]+= m18*u12+ m19*u13+ m20*u14+ m21*u15+ m22*u16+ m23*u17; wp[16]+= m24*u12+ m25*u13+ m26*u14+ m27*u15+ m28*u16+ m29*u17; wp[17]+= m30*u12+ m31*u13+ m32*u14+ m33*u15+ m34*u16+ m35*u17; wp[18]+= m0*u18 + m1*u19 + m2*u20 + m3*u21 + m4*u22 + m5*u23; wp[19]+= m6*u18 + m7*u19 + m8*u20 + m9*u21+ m10*u22+ m11*u23; wp[20]+= m12*u18+ m13*u19+ m14*u20+ m15*u21+ m16*u22+ m17*u23; wp[21]+= m18*u18+ m19*u19+ m20*u20+ m21*u21+ m22*u22+ m23*u23; wp[22]+= m24*u18+ m25*u19+ m26*u20+ m27*u21+ m28*u22+ m29*u23; wp[23]+= m30*u18+ m31*u19+ m32*u20+ m33*u21+ m34*u22+ m35*u23; wp[24]+= m0*u24 + m1*u25 + m2*u26 + m3*u27 + m4*u28 + m5*u29; wp[25]+= m6*u24 + m7*u25 + m8*u26 + m9*u27+ m10*u28+ m11*u29; wp[26]+= m12*u24+ m13*u25+ m14*u26+ m15*u27+ m16*u28+ m17*u29; wp[27]+= m18*u24+ m19*u25+ m20*u26+ m21*u27+ m22*u28+ m23*u29; wp[28]+= m24*u24+ m25*u25+ m26*u26+ m27*u27+ m28*u28+ m29*u29; wp[29]+= m30*u24+ m31*u25+ m32*u26+ m33*u27+ m34*u28+ m35*u29; wp[30]+= m0*u30 + m1*u31 + m2*u32 + m3*u33 + m4*u34 + m5*u35; wp[31]+= m6*u30 + m7*u31 + m8*u32 + m9*u33+ m10*u34+ m11*u35; wp[32]+= m12*u30+ m13*u31+ m14*u32+ m15*u33+ m16*u34+ m17*u35; wp[33]+= m18*u30+ m19*u31+ m20*u32+ m21*u33+ m22*u34+ m23*u35; wp[34]+= m24*u30+ m25*u31+ m26*u32+ m27*u33+ m28*u34+ m29*u35; wp[35]+= m30*u30+ m31*u31+ m32*u32+ m33*u33+ m34*u34+ m35*u35; } ierr = PetscLogFlops(2.0*216.0*(jmax-jmin));CHKERRQ(ierr); /* ... add i to row list for next nonzero entry */ il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; /* update jl */ } i = nexti; } /* save nonzero entries in k-th row of U ... */ /* invert diagonal block */ d = ba+k*36; ierr = PetscMemcpy(d,dk,36*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscKernel_A_gets_inverse_A_6(d,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; jmin = bi[k]; jmax = bi[k+1]; if (jmin < jmax) { for (j=jmin; j<jmax; j++) { vj = bj[j]; /* block col. index of U */ u = ba + j*36; wp = w + vj*36; for (k1=0; k1<36; k1++) { *u++ = *wp; *wp++ = 0.0; } } /* ... add k to row list for first nonzero entry in k-th row */ il[k] = jmin; i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; } } ierr = PetscFree(w);CHKERRQ(ierr); ierr = PetscFree2(il,jl);CHKERRQ(ierr); ierr = PetscFree2(dk,uik);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(1.3333*216*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); }