xref: /petsc/src/mat/impls/baij/seq/baijfact13.c (revision b2b2dd246975d7f9c8a1571def503d28e659d8b1)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
383287d42SBarry Smith /*
483287d42SBarry Smith     Factorization code for BAIJ format.
583287d42SBarry Smith */
67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
7c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
883287d42SBarry Smith 
983287d42SBarry Smith /*
1083287d42SBarry Smith       Version for when blocks are 3 by 3
1183287d42SBarry Smith */
124a2ae208SSatish Balay #undef __FUNCT__
134a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3"
140481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat C,Mat A,const MatFactorInfo *info)
1583287d42SBarry Smith {
1683287d42SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
1783287d42SBarry Smith   IS             isrow = b->row,isicol = b->icol;
186849ba73SBarry Smith   PetscErrorCode ierr;
195d0c19d7SBarry Smith   const PetscInt *r,*ic;
205d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
21690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j;
22690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,idx,*pj;
2383287d42SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
2483287d42SBarry Smith   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
2583287d42SBarry Smith   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
2683287d42SBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
2762bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
2883287d42SBarry Smith 
2983287d42SBarry Smith   PetscFunctionBegin;
3083287d42SBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
3183287d42SBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
32b0a32e0cSBarry Smith   ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
3383287d42SBarry Smith 
3483287d42SBarry Smith   for (i=0; i<n; i++) {
3583287d42SBarry Smith     nz    = bi[i+1] - bi[i];
3683287d42SBarry Smith     ajtmp = bj + bi[i];
3783287d42SBarry Smith     for  (j=0; j<nz; j++) {
3883287d42SBarry Smith       x = rtmp + 9*ajtmp[j];
3983287d42SBarry Smith       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0;
4083287d42SBarry Smith     }
4183287d42SBarry Smith     /* load in initial (unfactored row) */
4283287d42SBarry Smith     idx      = r[i];
4383287d42SBarry Smith     nz       = ai[idx+1] - ai[idx];
4483287d42SBarry Smith     ajtmpold = aj + ai[idx];
4583287d42SBarry Smith     v        = aa + 9*ai[idx];
4683287d42SBarry Smith     for (j=0; j<nz; j++) {
4783287d42SBarry Smith       x    = rtmp + 9*ic[ajtmpold[j]];
4883287d42SBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
4983287d42SBarry Smith       x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
5083287d42SBarry Smith       v    += 9;
5183287d42SBarry Smith     }
5283287d42SBarry Smith     row = *ajtmp++;
5383287d42SBarry Smith     while (row < i) {
5483287d42SBarry Smith       pc = rtmp + 9*row;
5583287d42SBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
5683287d42SBarry Smith       p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
5783287d42SBarry Smith       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
5883287d42SBarry Smith           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
5983287d42SBarry Smith         pv = ba + 9*diag_offset[row];
6083287d42SBarry Smith         pj = bj + diag_offset[row] + 1;
6183287d42SBarry Smith         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
6283287d42SBarry Smith         x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
6383287d42SBarry Smith         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
6483287d42SBarry Smith         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
6583287d42SBarry Smith         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
6683287d42SBarry Smith 
6783287d42SBarry Smith         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
6883287d42SBarry Smith         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
6983287d42SBarry Smith         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
7083287d42SBarry Smith 
7183287d42SBarry Smith         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
7283287d42SBarry Smith         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
7383287d42SBarry Smith         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
7483287d42SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
7583287d42SBarry Smith         pv += 9;
7683287d42SBarry Smith         for (j=0; j<nz; j++) {
7783287d42SBarry Smith           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
7883287d42SBarry Smith           x5   = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
7983287d42SBarry Smith           x    = rtmp + 9*pj[j];
8083287d42SBarry Smith           x[0] -= m1*x1 + m4*x2 + m7*x3;
8183287d42SBarry Smith           x[1] -= m2*x1 + m5*x2 + m8*x3;
8283287d42SBarry Smith           x[2] -= m3*x1 + m6*x2 + m9*x3;
8383287d42SBarry Smith 
8483287d42SBarry Smith           x[3] -= m1*x4 + m4*x5 + m7*x6;
8583287d42SBarry Smith           x[4] -= m2*x4 + m5*x5 + m8*x6;
8683287d42SBarry Smith           x[5] -= m3*x4 + m6*x5 + m9*x6;
8783287d42SBarry Smith 
8883287d42SBarry Smith           x[6] -= m1*x7 + m4*x8 + m7*x9;
8983287d42SBarry Smith           x[7] -= m2*x7 + m5*x8 + m8*x9;
9083287d42SBarry Smith           x[8] -= m3*x7 + m6*x8 + m9*x9;
9183287d42SBarry Smith           pv   += 9;
9283287d42SBarry Smith         }
93dc0b31edSSatish Balay         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
9483287d42SBarry Smith       }
9583287d42SBarry Smith       row = *ajtmp++;
9683287d42SBarry Smith     }
9783287d42SBarry Smith     /* finished row so stick it into b->a */
9883287d42SBarry Smith     pv = ba + 9*bi[i];
9983287d42SBarry Smith     pj = bj + bi[i];
10083287d42SBarry Smith     nz = bi[i+1] - bi[i];
10183287d42SBarry Smith     for (j=0; j<nz; j++) {
10283287d42SBarry Smith       x     = rtmp + 9*pj[j];
10383287d42SBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
10483287d42SBarry Smith       pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
10583287d42SBarry Smith       pv   += 9;
10683287d42SBarry Smith     }
10783287d42SBarry Smith     /* invert diagonal block */
10883287d42SBarry Smith     w = ba + 9*diag_offset[i];
10962bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr);
11083287d42SBarry Smith   }
11183287d42SBarry Smith 
11283287d42SBarry Smith   ierr = PetscFree(rtmp);CHKERRQ(ierr);
11383287d42SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
11483287d42SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
115db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqBAIJ_3;
116db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3;
11783287d42SBarry Smith   C->assembled = PETSC_TRUE;
118efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
11983287d42SBarry Smith   PetscFunctionReturn(0);
12083287d42SBarry Smith }
12117542729SShri Abhyankar 
12217542729SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_3_newdatastruct -
12317542729SShri Abhyankar      copied from MatLUFactorNumeric_SeqBAIJ_N_newdatastruct() and manually re-implemented
12417542729SShri Abhyankar        Kernel_A_gets_A_times_B()
12517542729SShri Abhyankar        Kernel_A_gets_A_minus_B_times_C()
12617542729SShri Abhyankar        Kernel_A_gets_inverse_A()
12717542729SShri Abhyankar */
12817542729SShri Abhyankar #undef __FUNCT__
12917542729SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_newdatastruct"
13017542729SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
13117542729SShri Abhyankar {
13217542729SShri Abhyankar   Mat            C=B;
13317542729SShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
13417542729SShri Abhyankar   IS             isrow = b->row,isicol = b->icol;
13517542729SShri Abhyankar   PetscErrorCode ierr;
13617542729SShri Abhyankar   const PetscInt *r,*ic,*ics;
13717542729SShri Abhyankar   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
13817542729SShri Abhyankar   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
13917542729SShri Abhyankar   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
14017542729SShri Abhyankar   PetscInt       bs2 = a->bs2,flg;
14117542729SShri Abhyankar   PetscReal      shift = info->shiftinblocks;
14217542729SShri Abhyankar 
14317542729SShri Abhyankar   PetscFunctionBegin;
14417542729SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
14517542729SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
14617542729SShri Abhyankar 
14717542729SShri Abhyankar   /* generate work space needed by the factorization */
14817542729SShri Abhyankar   ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
14917542729SShri Abhyankar   mwork = rtmp + bs2*n;
15017542729SShri Abhyankar   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
15117542729SShri Abhyankar   ics  = ic;
15217542729SShri Abhyankar 
15317542729SShri Abhyankar   for (i=0; i<n; i++){
15417542729SShri Abhyankar     /* zero rtmp */
15517542729SShri Abhyankar     /* L part */
15617542729SShri Abhyankar     nz    = bi[i+1] - bi[i];
15717542729SShri Abhyankar     bjtmp = bj + bi[i];
15817542729SShri Abhyankar     for  (j=0; j<nz; j++){
15917542729SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
16017542729SShri Abhyankar     }
16117542729SShri Abhyankar 
16217542729SShri Abhyankar     /* U part */
16317542729SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i];
16417542729SShri Abhyankar     bjtmp = bj + bi[2*n-i];
16517542729SShri Abhyankar     for  (j=0; j<nz; j++){
16617542729SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
16717542729SShri Abhyankar     }
16817542729SShri Abhyankar 
16917542729SShri Abhyankar     /* load in initial (unfactored row) */
17017542729SShri Abhyankar     nz    = ai[r[i]+1] - ai[r[i]];
17117542729SShri Abhyankar     ajtmp = aj + ai[r[i]];
17217542729SShri Abhyankar     v     = aa + bs2*ai[r[i]];
17317542729SShri Abhyankar     for (j=0; j<nz; j++) {
17417542729SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
17517542729SShri Abhyankar     }
17617542729SShri Abhyankar 
17717542729SShri Abhyankar     /* elimination */
17817542729SShri Abhyankar     bjtmp = bj + bi[i];
17917542729SShri Abhyankar     nzL   = bi[i+1] - bi[i];
180b1646270SShri Abhyankar     for(k = 0;k < nzL;k++){
181b1646270SShri Abhyankar       row = bjtmp[k];
18217542729SShri Abhyankar       pc = rtmp + bs2*row;
18317542729SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
18417542729SShri Abhyankar       if (flg) {
18517542729SShri Abhyankar         pv = b->a + bs2*bdiag[row];
18617542729SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
18717542729SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
18817542729SShri Abhyankar 
18917542729SShri Abhyankar         pj = b->j + bi[2*n-row]; /* begining of U(row,:) */
19017542729SShri Abhyankar         pv = b->a + bs2*bi[2*n-row];
19117542729SShri Abhyankar         nz = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */
19217542729SShri Abhyankar         for (j=0; j<nz; j++) {
19317542729SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
19417542729SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
19517542729SShri Abhyankar           v    = rtmp + bs2*pj[j];
19617542729SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
19717542729SShri Abhyankar           pv  += bs2;
19817542729SShri Abhyankar         }
19917542729SShri Abhyankar         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
20017542729SShri Abhyankar       }
20117542729SShri Abhyankar     }
20217542729SShri Abhyankar 
20317542729SShri Abhyankar     /* finished row so stick it into b->a */
20417542729SShri Abhyankar     /* L part */
20517542729SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
20617542729SShri Abhyankar     pj   = b->j + bi[i] ;
20717542729SShri Abhyankar     nz   = bi[i+1] - bi[i];
20817542729SShri Abhyankar     for (j=0; j<nz; j++) {
20917542729SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
21017542729SShri Abhyankar     }
21117542729SShri Abhyankar 
21217542729SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
21317542729SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
21417542729SShri Abhyankar     pj   = b->j + bdiag[i];
21517542729SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
21617542729SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
21717542729SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr);
21817542729SShri Abhyankar 
21917542729SShri Abhyankar     /* U part */
22017542729SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
22117542729SShri Abhyankar     pj = b->j + bi[2*n-i];
22217542729SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
22317542729SShri Abhyankar     for (j=0; j<nz; j++){
22417542729SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
22517542729SShri Abhyankar     }
22617542729SShri Abhyankar   }
22717542729SShri Abhyankar 
22817542729SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
22917542729SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
23017542729SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
23117542729SShri Abhyankar 
23217542729SShri Abhyankar   C->assembled = PETSC_TRUE;
23317542729SShri Abhyankar   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
23417542729SShri Abhyankar   PetscFunctionReturn(0);
23517542729SShri Abhyankar }
23617542729SShri Abhyankar 
23717542729SShri Abhyankar #undef __FUNCT__
23817542729SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering"
23917542729SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
24017542729SShri Abhyankar {
24117542729SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
24217542729SShri Abhyankar   PetscErrorCode ierr;
24317542729SShri Abhyankar   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
24417542729SShri Abhyankar   PetscInt       *ajtmpold,*ajtmp,nz,row;
24517542729SShri Abhyankar   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
24617542729SShri Abhyankar   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
24717542729SShri Abhyankar   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
24817542729SShri Abhyankar   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
24917542729SShri Abhyankar   MatScalar      *ba = b->a,*aa = a->a;
25017542729SShri Abhyankar   PetscReal      shift = info->shiftinblocks;
25117542729SShri Abhyankar 
25217542729SShri Abhyankar   PetscFunctionBegin;
25317542729SShri Abhyankar   ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
25417542729SShri Abhyankar 
25517542729SShri Abhyankar   for (i=0; i<n; i++) {
25617542729SShri Abhyankar     nz    = bi[i+1] - bi[i];
25717542729SShri Abhyankar     ajtmp = bj + bi[i];
25817542729SShri Abhyankar     for  (j=0; j<nz; j++) {
25917542729SShri Abhyankar       x = rtmp+9*ajtmp[j];
26017542729SShri Abhyankar       x[0]  = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = 0.0;
26117542729SShri Abhyankar     }
26217542729SShri Abhyankar     /* load in initial (unfactored row) */
26317542729SShri Abhyankar     nz       = ai[i+1] - ai[i];
26417542729SShri Abhyankar     ajtmpold = aj + ai[i];
26517542729SShri Abhyankar     v        = aa + 9*ai[i];
26617542729SShri Abhyankar     for (j=0; j<nz; j++) {
26717542729SShri Abhyankar       x    = rtmp+9*ajtmpold[j];
26817542729SShri Abhyankar       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
26917542729SShri Abhyankar       x[4]  = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
27017542729SShri Abhyankar       v    += 9;
27117542729SShri Abhyankar     }
27217542729SShri Abhyankar     row = *ajtmp++;
27317542729SShri Abhyankar     while (row < i) {
27417542729SShri Abhyankar       pc  = rtmp + 9*row;
27517542729SShri Abhyankar       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
27617542729SShri Abhyankar       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
27717542729SShri Abhyankar       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
27817542729SShri Abhyankar           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
27917542729SShri Abhyankar         pv = ba + 9*diag_offset[row];
28017542729SShri Abhyankar         pj = bj + diag_offset[row] + 1;
28117542729SShri Abhyankar         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
28217542729SShri Abhyankar         x5  = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
28317542729SShri Abhyankar         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
28417542729SShri Abhyankar         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
28517542729SShri Abhyankar         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
28617542729SShri Abhyankar 
28717542729SShri Abhyankar         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
28817542729SShri Abhyankar         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
28917542729SShri Abhyankar         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
29017542729SShri Abhyankar 
29117542729SShri Abhyankar         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
29217542729SShri Abhyankar         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
29317542729SShri Abhyankar         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
29417542729SShri Abhyankar 
29517542729SShri Abhyankar         nz = bi[row+1] - diag_offset[row] - 1;
29617542729SShri Abhyankar         pv += 9;
29717542729SShri Abhyankar         for (j=0; j<nz; j++) {
29817542729SShri Abhyankar           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
29917542729SShri Abhyankar           x5   = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
30017542729SShri Abhyankar           x    = rtmp + 9*pj[j];
30117542729SShri Abhyankar           x[0] -= m1*x1 + m4*x2 + m7*x3;
30217542729SShri Abhyankar           x[1] -= m2*x1 + m5*x2 + m8*x3;
30317542729SShri Abhyankar           x[2] -= m3*x1 + m6*x2 + m9*x3;
30417542729SShri Abhyankar 
30517542729SShri Abhyankar           x[3] -= m1*x4 + m4*x5 + m7*x6;
30617542729SShri Abhyankar           x[4] -= m2*x4 + m5*x5 + m8*x6;
30717542729SShri Abhyankar           x[5] -= m3*x4 + m6*x5 + m9*x6;
30817542729SShri Abhyankar 
30917542729SShri Abhyankar           x[6] -= m1*x7 + m4*x8 + m7*x9;
31017542729SShri Abhyankar           x[7] -= m2*x7 + m5*x8 + m8*x9;
31117542729SShri Abhyankar           x[8] -= m3*x7 + m6*x8 + m9*x9;
31217542729SShri Abhyankar           pv   += 9;
31317542729SShri Abhyankar         }
31417542729SShri Abhyankar         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
31517542729SShri Abhyankar       }
31617542729SShri Abhyankar       row = *ajtmp++;
31717542729SShri Abhyankar     }
31817542729SShri Abhyankar     /* finished row so stick it into b->a */
31917542729SShri Abhyankar     pv = ba + 9*bi[i];
32017542729SShri Abhyankar     pj = bj + bi[i];
32117542729SShri Abhyankar     nz = bi[i+1] - bi[i];
32217542729SShri Abhyankar     for (j=0; j<nz; j++) {
32317542729SShri Abhyankar       x      = rtmp+9*pj[j];
32417542729SShri Abhyankar       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
32517542729SShri Abhyankar       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
32617542729SShri Abhyankar       pv   += 9;
32717542729SShri Abhyankar     }
32817542729SShri Abhyankar     /* invert diagonal block */
32917542729SShri Abhyankar     w = ba + 9*diag_offset[i];
33017542729SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr);
33117542729SShri Abhyankar   }
33217542729SShri Abhyankar 
33317542729SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
33417542729SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_3_NaturalOrdering;
33517542729SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
33617542729SShri Abhyankar   C->assembled = PETSC_TRUE;
33717542729SShri Abhyankar   ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
33817542729SShri Abhyankar   PetscFunctionReturn(0);
33917542729SShri Abhyankar }
34017542729SShri Abhyankar 
34117542729SShri Abhyankar /*
34217542729SShri Abhyankar   MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct -
34317542729SShri Abhyankar     copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct()
34417542729SShri Abhyankar */
34517542729SShri Abhyankar #undef __FUNCT__
34617542729SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct"
34717542729SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
34817542729SShri Abhyankar {
34917542729SShri Abhyankar   Mat            C=B;
35017542729SShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
35117542729SShri Abhyankar   PetscErrorCode ierr;
35217542729SShri Abhyankar   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
35317542729SShri Abhyankar   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
35417542729SShri Abhyankar   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
35517542729SShri Abhyankar   PetscInt       bs2 = a->bs2,flg;
35617542729SShri Abhyankar   PetscReal      shift = info->shiftinblocks;
35717542729SShri Abhyankar 
35817542729SShri Abhyankar   PetscFunctionBegin;
35917542729SShri Abhyankar   /* generate work space needed by the factorization */
36017542729SShri Abhyankar   ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
36117542729SShri Abhyankar   mwork = rtmp + bs2*n;
36217542729SShri Abhyankar   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
36317542729SShri Abhyankar 
36417542729SShri Abhyankar   for (i=0; i<n; i++){
36517542729SShri Abhyankar     /* zero rtmp */
36617542729SShri Abhyankar     /* L part */
36717542729SShri Abhyankar     nz    = bi[i+1] - bi[i];
36817542729SShri Abhyankar     bjtmp = bj + bi[i];
36917542729SShri Abhyankar     for  (j=0; j<nz; j++){
37017542729SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
37117542729SShri Abhyankar     }
37217542729SShri Abhyankar 
37317542729SShri Abhyankar     /* U part */
37417542729SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i];
37517542729SShri Abhyankar     bjtmp = bj + bi[2*n-i];
37617542729SShri Abhyankar     for  (j=0; j<nz; j++){
37717542729SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
37817542729SShri Abhyankar     }
37917542729SShri Abhyankar 
38017542729SShri Abhyankar     /* load in initial (unfactored row) */
38117542729SShri Abhyankar     nz    = ai[i+1] - ai[i];
38217542729SShri Abhyankar     ajtmp = aj + ai[i];
38317542729SShri Abhyankar     v     = aa + bs2*ai[i];
38417542729SShri Abhyankar     for (j=0; j<nz; j++) {
38517542729SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
38617542729SShri Abhyankar     }
38717542729SShri Abhyankar 
38817542729SShri Abhyankar     /* elimination */
38917542729SShri Abhyankar     bjtmp = bj + bi[i];
39017542729SShri Abhyankar     nzL   = bi[i+1] - bi[i];
391b1646270SShri Abhyankar     for(k=0;k<nzL;k++){
392b1646270SShri Abhyankar       row = bjtmp[k];
39317542729SShri Abhyankar       pc = rtmp + bs2*row;
39417542729SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
39517542729SShri Abhyankar       if (flg) {
39617542729SShri Abhyankar         pv = b->a + bs2*bdiag[row];
39717542729SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
39817542729SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
39917542729SShri Abhyankar 
40017542729SShri Abhyankar         pj = b->j + bi[2*n-row]; /* begining of U(row,:) */
40117542729SShri Abhyankar         pv = b->a + bs2*bi[2*n-row];
40217542729SShri Abhyankar         nz = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */
40317542729SShri Abhyankar         for (j=0; j<nz; j++) {
40417542729SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
40517542729SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
40617542729SShri Abhyankar           v    = rtmp + bs2*pj[j];
40717542729SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
40817542729SShri Abhyankar           pv  += bs2;
40917542729SShri Abhyankar         }
41017542729SShri Abhyankar         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
41117542729SShri Abhyankar       }
41217542729SShri Abhyankar     }
41317542729SShri Abhyankar 
41417542729SShri Abhyankar     /* finished row so stick it into b->a */
41517542729SShri Abhyankar     /* L part */
41617542729SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
41717542729SShri Abhyankar     pj   = b->j + bi[i] ;
41817542729SShri Abhyankar     nz   = bi[i+1] - bi[i];
41917542729SShri Abhyankar     for (j=0; j<nz; j++) {
42017542729SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
42117542729SShri Abhyankar     }
42217542729SShri Abhyankar 
42317542729SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
42417542729SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
42517542729SShri Abhyankar     pj   = b->j + bdiag[i];
42617542729SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
42717542729SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
42817542729SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr);
42917542729SShri Abhyankar 
43017542729SShri Abhyankar     /* U part */
43117542729SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
43217542729SShri Abhyankar     pj = b->j + bi[2*n-i];
43317542729SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
43417542729SShri Abhyankar     for (j=0; j<nz; j++){
43517542729SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
43617542729SShri Abhyankar     }
43717542729SShri Abhyankar   }
43817542729SShri Abhyankar 
43917542729SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
44017542729SShri Abhyankar   C->assembled = PETSC_TRUE;
44117542729SShri Abhyankar   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
44217542729SShri Abhyankar   PetscFunctionReturn(0);
44317542729SShri Abhyankar }
444*b2b2dd24SShri Abhyankar 
445*b2b2dd24SShri Abhyankar #undef __FUNCT__
446*b2b2dd24SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct_v2"
447*b2b2dd24SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct_v2(Mat B,Mat A,const MatFactorInfo *info)
448*b2b2dd24SShri Abhyankar {
449*b2b2dd24SShri Abhyankar   Mat            C=B;
450*b2b2dd24SShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
451*b2b2dd24SShri Abhyankar   PetscErrorCode ierr;
452*b2b2dd24SShri Abhyankar   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
453*b2b2dd24SShri Abhyankar   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
454*b2b2dd24SShri Abhyankar   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
455*b2b2dd24SShri Abhyankar   PetscInt       bs2 = a->bs2,flg;
456*b2b2dd24SShri Abhyankar   PetscReal      shift = info->shiftinblocks;
457*b2b2dd24SShri Abhyankar 
458*b2b2dd24SShri Abhyankar   PetscFunctionBegin;
459*b2b2dd24SShri Abhyankar   /* generate work space needed by the factorization */
460*b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
461*b2b2dd24SShri Abhyankar   mwork = rtmp + bs2*n;
462*b2b2dd24SShri Abhyankar   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
463*b2b2dd24SShri Abhyankar 
464*b2b2dd24SShri Abhyankar   for (i=0; i<n; i++){
465*b2b2dd24SShri Abhyankar     /* zero rtmp */
466*b2b2dd24SShri Abhyankar     /* L part */
467*b2b2dd24SShri Abhyankar     nz    = bi[i+1] - bi[i];
468*b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
469*b2b2dd24SShri Abhyankar     for  (j=0; j<nz; j++){
470*b2b2dd24SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
471*b2b2dd24SShri Abhyankar     }
472*b2b2dd24SShri Abhyankar 
473*b2b2dd24SShri Abhyankar     /* U part */
474*b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1];
475*b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i+1] + 1;
476*b2b2dd24SShri Abhyankar     for  (j=0; j<nz; j++){
477*b2b2dd24SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
478*b2b2dd24SShri Abhyankar     }
479*b2b2dd24SShri Abhyankar 
480*b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
481*b2b2dd24SShri Abhyankar     nz    = ai[i+1] - ai[i];
482*b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
483*b2b2dd24SShri Abhyankar     v     = aa + bs2*ai[i];
484*b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
485*b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
486*b2b2dd24SShri Abhyankar     }
487*b2b2dd24SShri Abhyankar 
488*b2b2dd24SShri Abhyankar     /* elimination */
489*b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
490*b2b2dd24SShri Abhyankar     nzL   = bi[i+1] - bi[i];
491*b2b2dd24SShri Abhyankar     for(k=0;k<nzL;k++){
492*b2b2dd24SShri Abhyankar       row = bjtmp[k];
493*b2b2dd24SShri Abhyankar       pc = rtmp + bs2*row;
494*b2b2dd24SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
495*b2b2dd24SShri Abhyankar       if (flg) {
496*b2b2dd24SShri Abhyankar         pv = b->a + bs2*bdiag[row];
497*b2b2dd24SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
498*b2b2dd24SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
499*b2b2dd24SShri Abhyankar 
500*b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
501*b2b2dd24SShri Abhyankar 	pv = b->a + bs2*(bdiag[row+1]+1);
502*b2b2dd24SShri Abhyankar 	nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
503*b2b2dd24SShri Abhyankar         for (j=0; j<nz; j++) {
504*b2b2dd24SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
505*b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
506*b2b2dd24SShri Abhyankar           v    = rtmp + bs2*pj[j];
507*b2b2dd24SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
508*b2b2dd24SShri Abhyankar           pv  += bs2;
509*b2b2dd24SShri Abhyankar         }
510*b2b2dd24SShri Abhyankar         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
511*b2b2dd24SShri Abhyankar       }
512*b2b2dd24SShri Abhyankar     }
513*b2b2dd24SShri Abhyankar 
514*b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
515*b2b2dd24SShri Abhyankar     /* L part */
516*b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
517*b2b2dd24SShri Abhyankar     pj   = b->j + bi[i] ;
518*b2b2dd24SShri Abhyankar     nz   = bi[i+1] - bi[i];
519*b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
520*b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
521*b2b2dd24SShri Abhyankar     }
522*b2b2dd24SShri Abhyankar 
523*b2b2dd24SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
524*b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
525*b2b2dd24SShri Abhyankar     pj   = b->j + bdiag[i];
526*b2b2dd24SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
527*b2b2dd24SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
528*b2b2dd24SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr);
529*b2b2dd24SShri Abhyankar 
530*b2b2dd24SShri Abhyankar     /* U part */
531*b2b2dd24SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
532*b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
533*b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
534*b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++){
535*b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
536*b2b2dd24SShri Abhyankar     }
537*b2b2dd24SShri Abhyankar   }
538*b2b2dd24SShri Abhyankar 
539*b2b2dd24SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
540*b2b2dd24SShri Abhyankar   C->assembled = PETSC_TRUE;
541*b2b2dd24SShri Abhyankar   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
542*b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
543*b2b2dd24SShri Abhyankar }
544*b2b2dd24SShri Abhyankar 
545