xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact9.c (revision 719d5645761d844e4357b7ee00a3296dffe0b787)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
33a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
481278733SSatish Balay #include "src/inline/ilu.h"
581278733SSatish Balay 
681278733SSatish Balay /* Version for when blocks are 6 by 6 */
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_6"
9*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat C,Mat A,MatFactorInfo *info)
1081278733SSatish Balay {
1181278733SSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
1281278733SSatish Balay   IS             perm = b->row;
136849ba73SBarry Smith   PetscErrorCode ierr;
1413f74950SBarry Smith   PetscInt       *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
1513f74950SBarry Smith   PetscInt       *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1681278733SSatish Balay   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
171b3064deSBarry Smith   MatScalar      *u,*d,*w,*wp,u0,u1,u2,u3,u4,u5,u6,u7,u8,u9,u10,u11,u12;
181b3064deSBarry Smith   MatScalar      u13,u14,u15,u16,u17,u18,u19,u20,u21,u22,u23,u24,u25,u26,u27;
191b3064deSBarry Smith   MatScalar      u28,u29,u30,u31,u32,u33,u34,u35;
2062bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
2181278733SSatish Balay 
2281278733SSatish Balay   PetscFunctionBegin;
2381278733SSatish Balay   /* initialization */
2481278733SSatish Balay   ierr = PetscMalloc(36*mbs*sizeof(MatScalar),&w);CHKERRQ(ierr);
2581278733SSatish Balay   ierr = PetscMemzero(w,36*mbs*sizeof(MatScalar));CHKERRQ(ierr);
2613f74950SBarry Smith   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
2781278733SSatish Balay   jl = il + mbs;
2881278733SSatish Balay   for (i=0; i<mbs; i++) {
2981278733SSatish Balay     jl[i] = mbs; il[0] = 0;
3081278733SSatish Balay   }
3181278733SSatish Balay   ierr = PetscMalloc(72*sizeof(MatScalar),&dk);CHKERRQ(ierr);
3281278733SSatish Balay   uik  = dk + 36;
3381278733SSatish Balay   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
3481278733SSatish Balay 
3581278733SSatish Balay   /* check permutation */
3681278733SSatish Balay   if (!a->permute){
3781278733SSatish Balay     ai = a->i; aj = a->j; aa = a->a;
3881278733SSatish Balay   } else {
3981278733SSatish Balay     ai   = a->inew; aj = a->jnew;
4081278733SSatish Balay     ierr = PetscMalloc(36*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
4181278733SSatish Balay     ierr = PetscMemcpy(aa,a->a,36*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
4213f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
4313f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
4481278733SSatish Balay 
4581278733SSatish Balay     for (i=0; i<mbs; i++){
4681278733SSatish Balay       jmin = ai[i]; jmax = ai[i+1];
4781278733SSatish Balay       for (j=jmin; j<jmax; j++){
4881278733SSatish Balay         while (a2anew[j] != j){
4981278733SSatish Balay           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
5081278733SSatish Balay           for (k1=0; k1<36; k1++){
5181278733SSatish Balay             dk[k1]       = aa[k*36+k1];
5281278733SSatish Balay             aa[k*36+k1] = aa[j*36+k1];
5381278733SSatish Balay             aa[j*36+k1] = dk[k1];
5481278733SSatish Balay           }
5581278733SSatish Balay         }
5681278733SSatish Balay         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
5781278733SSatish Balay         if (i > aj[j]){
5881278733SSatish Balay           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
5981278733SSatish Balay           ap = aa + j*36;                     /* ptr to the beginning of j-th block of aa */
6081278733SSatish Balay           for (k=0; k<36; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6181278733SSatish Balay           for (k=0; k<6; k++){               /* j-th block of aa <- dk^T */
6281278733SSatish Balay             for (k1=0; k1<6; k1++) *ap++ = dk[k + 6*k1];
6381278733SSatish Balay           }
6481278733SSatish Balay         }
6581278733SSatish Balay       }
6681278733SSatish Balay     }
6781278733SSatish Balay     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6881278733SSatish Balay   }
6981278733SSatish Balay 
7081278733SSatish Balay   /* for each row k */
7181278733SSatish Balay   for (k = 0; k<mbs; k++){
7281278733SSatish Balay 
7381278733SSatish Balay     /*initialize k-th row with elements nonzero in row perm(k) of A */
7481278733SSatish Balay     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
7581278733SSatish Balay     if (jmin < jmax) {
7681278733SSatish Balay       ap = aa + jmin*36;
7781278733SSatish Balay       for (j = jmin; j < jmax; j++){
7881278733SSatish Balay         vj = perm_ptr[aj[j]];         /* block col. index */
7981278733SSatish Balay         wp = w + vj*36;
8081278733SSatish Balay         for (i=0; i<36; i++) *wp++ = *ap++;
8181278733SSatish Balay       }
8281278733SSatish Balay     }
8381278733SSatish Balay 
8481278733SSatish Balay     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
8581278733SSatish Balay     ierr = PetscMemcpy(dk,w+k*36,36*sizeof(MatScalar));CHKERRQ(ierr);
8681278733SSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
8781278733SSatish Balay 
8881278733SSatish Balay     while (i < mbs){
8981278733SSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
9081278733SSatish Balay 
9181278733SSatish Balay       /* compute multiplier */
9281278733SSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9381278733SSatish Balay 
9481278733SSatish Balay       /* uik = -inv(Di)*U_bar(i,k) */
9581278733SSatish Balay       d = ba + i*36;
9681278733SSatish Balay       u = ba + ili*36;
9781278733SSatish Balay 
981b3064deSBarry Smith       u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6];
991b3064deSBarry Smith       u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13];
1001b3064deSBarry Smith       u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20];
1011b3064deSBarry Smith       u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27];
1021b3064deSBarry Smith       u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34];
1031b3064deSBarry Smith       u35 = u[35];
10481278733SSatish Balay 
1051b3064deSBarry Smith       uik[0] = -(d[0]*u0 + d[6]*u1 + d[12]*u2 + d[18]*u3 + d[24]*u4 + d[30]*u5);
1061b3064deSBarry Smith       uik[1] = -(d[1]*u0 + d[7]*u1 + d[13]*u2 + d[19]*u3 + d[25]*u4 + d[31]*u5);
1071b3064deSBarry Smith       uik[2] = -(d[2]*u0 + d[8]*u1 + d[14]*u2 + d[20]*u3 + d[26]*u4 + d[32]*u5);
1081b3064deSBarry Smith       uik[3] = -(d[3]*u0 + d[9]*u1 + d[15]*u2 + d[21]*u3 + d[27]*u4 + d[33]*u5);
1091b3064deSBarry Smith       uik[4] = -(d[4]*u0+ d[10]*u1 + d[16]*u2 + d[22]*u3 + d[28]*u4 + d[34]*u5);
1101b3064deSBarry Smith       uik[5] = -(d[5]*u0+ d[11]*u1 + d[17]*u2 + d[23]*u3 + d[29]*u4 + d[35]*u5);
11181278733SSatish Balay 
1121b3064deSBarry Smith       uik[6] = -(d[0]*u6 + d[6]*u7 + d[12]*u8 + d[18]*u9 + d[24]*u10 + d[30]*u11);
1131b3064deSBarry Smith       uik[7] = -(d[1]*u6 + d[7]*u7 + d[13]*u8 + d[19]*u9 + d[25]*u10 + d[31]*u11);
1141b3064deSBarry Smith       uik[8] = -(d[2]*u6 + d[8]*u7 + d[14]*u8 + d[20]*u9 + d[26]*u10 + d[32]*u11);
1151b3064deSBarry Smith       uik[9] = -(d[3]*u6 + d[9]*u7 + d[15]*u8 + d[21]*u9 + d[27]*u10 + d[33]*u11);
1161b3064deSBarry Smith       uik[10]= -(d[4]*u6+ d[10]*u7 + d[16]*u8 + d[22]*u9 + d[28]*u10 + d[34]*u11);
1171b3064deSBarry Smith       uik[11]= -(d[5]*u6+ d[11]*u7 + d[17]*u8 + d[23]*u9 + d[29]*u10 + d[35]*u11);
11881278733SSatish Balay 
1191b3064deSBarry Smith       uik[12] = -(d[0]*u12 + d[6]*u13 + d[12]*u14 + d[18]*u15 + d[24]*u16 + d[30]*u17);
1201b3064deSBarry Smith       uik[13] = -(d[1]*u12 + d[7]*u13 + d[13]*u14 + d[19]*u15 + d[25]*u16 + d[31]*u17);
1211b3064deSBarry Smith       uik[14] = -(d[2]*u12 + d[8]*u13 + d[14]*u14 + d[20]*u15 + d[26]*u16 + d[32]*u17);
1221b3064deSBarry Smith       uik[15] = -(d[3]*u12 + d[9]*u13 + d[15]*u14 + d[21]*u15 + d[27]*u16 + d[33]*u17);
1231b3064deSBarry Smith       uik[16] = -(d[4]*u12+ d[10]*u13 + d[16]*u14 + d[22]*u15 + d[28]*u16 + d[34]*u17);
1241b3064deSBarry Smith       uik[17] = -(d[5]*u12+ d[11]*u13 + d[17]*u14 + d[23]*u15 + d[29]*u16 + d[35]*u17);
12581278733SSatish Balay 
1261b3064deSBarry Smith       uik[18] = -(d[0]*u18 + d[6]*u19 + d[12]*u20 + d[18]*u21 + d[24]*u22 + d[30]*u23);
1271b3064deSBarry Smith       uik[19] = -(d[1]*u18 + d[7]*u19 + d[13]*u20 + d[19]*u21 + d[25]*u22 + d[31]*u23);
1281b3064deSBarry Smith       uik[20] = -(d[2]*u18 + d[8]*u19 + d[14]*u20 + d[20]*u21 + d[26]*u22 + d[32]*u23);
1291b3064deSBarry Smith       uik[21] = -(d[3]*u18 + d[9]*u19 + d[15]*u20 + d[21]*u21 + d[27]*u22 + d[33]*u23);
1301b3064deSBarry Smith       uik[22] = -(d[4]*u18+ d[10]*u19 + d[16]*u20 + d[22]*u21 + d[28]*u22 + d[34]*u23);
1311b3064deSBarry Smith       uik[23] = -(d[5]*u18+ d[11]*u19 + d[17]*u20 + d[23]*u21 + d[29]*u22 + d[35]*u23);
13281278733SSatish Balay 
1331b3064deSBarry Smith       uik[24] = -(d[0]*u24 + d[6]*u25 + d[12]*u26 + d[18]*u27 + d[24]*u28 + d[30]*u29);
1341b3064deSBarry Smith       uik[25] = -(d[1]*u24 + d[7]*u25 + d[13]*u26 + d[19]*u27 + d[25]*u28 + d[31]*u29);
1351b3064deSBarry Smith       uik[26] = -(d[2]*u24 + d[8]*u25 + d[14]*u26 + d[20]*u27 + d[26]*u28 + d[32]*u29);
1361b3064deSBarry Smith       uik[27] = -(d[3]*u24 + d[9]*u25 + d[15]*u26 + d[21]*u27 + d[27]*u28 + d[33]*u29);
1371b3064deSBarry Smith       uik[28] = -(d[4]*u24+ d[10]*u25 + d[16]*u26 + d[22]*u27 + d[28]*u28 + d[34]*u29);
1381b3064deSBarry Smith       uik[29] = -(d[5]*u24+ d[11]*u25 + d[17]*u26 + d[23]*u27 + d[29]*u28 + d[35]*u29);
1391b3064deSBarry Smith 
1401b3064deSBarry Smith       uik[30] = -(d[0]*u30 + d[6]*u31 + d[12]*u32 + d[18]*u33 + d[24]*u34 + d[30]*u35);
1411b3064deSBarry Smith       uik[31] = -(d[1]*u30 + d[7]*u31 + d[13]*u32 + d[19]*u33 + d[25]*u34 + d[31]*u35);
1421b3064deSBarry Smith       uik[32] = -(d[2]*u30 + d[8]*u31 + d[14]*u32 + d[20]*u33 + d[26]*u34 + d[32]*u35);
1431b3064deSBarry Smith       uik[33] = -(d[3]*u30 + d[9]*u31 + d[15]*u32 + d[21]*u33 + d[27]*u34 + d[33]*u35);
1441b3064deSBarry Smith       uik[34] = -(d[4]*u30+ d[10]*u31 + d[16]*u32 + d[22]*u33 + d[28]*u34 + d[34]*u35);
1451b3064deSBarry Smith       uik[35] = -(d[5]*u30+ d[11]*u31 + d[17]*u32 + d[23]*u33 + d[29]*u34 + d[35]*u35);
14681278733SSatish Balay 
14781278733SSatish Balay       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
1481b3064deSBarry Smith       dk[0] +=  uik[0]*u0 + uik[1]*u1 + uik[2]*u2 + uik[3]*u3 + uik[4]*u4 + uik[5]*u5;
1491b3064deSBarry Smith       dk[1] +=  uik[6]*u0 + uik[7]*u1 + uik[8]*u2 + uik[9]*u3+ uik[10]*u4+ uik[11]*u5;
1501b3064deSBarry Smith       dk[2] += uik[12]*u0+ uik[13]*u1+ uik[14]*u2+ uik[15]*u3+ uik[16]*u4+ uik[17]*u5;
1511b3064deSBarry Smith       dk[3] += uik[18]*u0+ uik[19]*u1+ uik[20]*u2+ uik[21]*u3+ uik[22]*u4+ uik[23]*u5;
1521b3064deSBarry Smith       dk[4] += uik[24]*u0+ uik[25]*u1+ uik[26]*u2+ uik[27]*u3+ uik[28]*u4+ uik[29]*u5;
1531b3064deSBarry Smith       dk[5] += uik[30]*u0+ uik[31]*u1+ uik[32]*u2+ uik[33]*u3+ uik[34]*u4+ uik[35]*u5;
15481278733SSatish Balay 
1551b3064deSBarry Smith       dk[6] +=  uik[0]*u6 + uik[1]*u7 + uik[2]*u8 + uik[3]*u9 + uik[4]*u10 + uik[5]*u11;
1561b3064deSBarry Smith       dk[7] +=  uik[6]*u6 + uik[7]*u7 + uik[8]*u8 + uik[9]*u9+ uik[10]*u10+ uik[11]*u11;
1571b3064deSBarry Smith       dk[8] += uik[12]*u6+ uik[13]*u7+ uik[14]*u8+ uik[15]*u9+ uik[16]*u10+ uik[17]*u11;
1581b3064deSBarry Smith       dk[9] += uik[18]*u6+ uik[19]*u7+ uik[20]*u8+ uik[21]*u9+ uik[22]*u10+ uik[23]*u11;
1591b3064deSBarry Smith       dk[10]+= uik[24]*u6+ uik[25]*u7+ uik[26]*u8+ uik[27]*u9+ uik[28]*u10+ uik[29]*u11;
1601b3064deSBarry Smith       dk[11]+= uik[30]*u6+ uik[31]*u7+ uik[32]*u8+ uik[33]*u9+ uik[34]*u10+ uik[35]*u11;
16181278733SSatish Balay 
1621b3064deSBarry Smith       dk[12]+=  uik[0]*u12 + uik[1]*u13 + uik[2]*u14 + uik[3]*u15 + uik[4]*u16 + uik[5]*u17;
1631b3064deSBarry Smith       dk[13]+=  uik[6]*u12 + uik[7]*u13 + uik[8]*u14 + uik[9]*u15+ uik[10]*u16+ uik[11]*u17;
1641b3064deSBarry Smith       dk[14]+= uik[12]*u12+ uik[13]*u13+ uik[14]*u14+ uik[15]*u15+ uik[16]*u16+ uik[17]*u17;
1651b3064deSBarry Smith       dk[15]+= uik[18]*u12+ uik[19]*u13+ uik[20]*u14+ uik[21]*u15+ uik[22]*u16+ uik[23]*u17;
1661b3064deSBarry Smith       dk[16]+= uik[24]*u12+ uik[25]*u13+ uik[26]*u14+ uik[27]*u15+ uik[28]*u16+ uik[29]*u17;
1671b3064deSBarry Smith       dk[17]+= uik[30]*u12+ uik[31]*u13+ uik[32]*u14+ uik[33]*u15+ uik[34]*u16+ uik[35]*u17;
16881278733SSatish Balay 
1691b3064deSBarry Smith       dk[18]+=  uik[0]*u18 + uik[1]*u19 + uik[2]*u20 + uik[3]*u21 + uik[4]*u22 + uik[5]*u23;
1701b3064deSBarry Smith       dk[19]+=  uik[6]*u18 + uik[7]*u19 + uik[8]*u20 + uik[9]*u21+ uik[10]*u22+ uik[11]*u23;
1711b3064deSBarry Smith       dk[20]+= uik[12]*u18+ uik[13]*u19+ uik[14]*u20+ uik[15]*u21+ uik[16]*u22+ uik[17]*u23;
1721b3064deSBarry Smith       dk[21]+= uik[18]*u18+ uik[19]*u19+ uik[20]*u20+ uik[21]*u21+ uik[22]*u22+ uik[23]*u23;
1731b3064deSBarry Smith       dk[22]+= uik[24]*u18+ uik[25]*u19+ uik[26]*u20+ uik[27]*u21+ uik[28]*u22+ uik[29]*u23;
1741b3064deSBarry Smith       dk[23]+= uik[30]*u18+ uik[31]*u19+ uik[32]*u20+ uik[33]*u21+ uik[34]*u22+ uik[35]*u23;
17581278733SSatish Balay 
1761b3064deSBarry Smith       dk[24]+=  uik[0]*u24 + uik[1]*u25 + uik[2]*u26 + uik[3]*u27 + uik[4]*u28 + uik[5]*u29;
1771b3064deSBarry Smith       dk[25]+=  uik[6]*u24 + uik[7]*u25 + uik[8]*u26 + uik[9]*u27+ uik[10]*u28+ uik[11]*u29;
1781b3064deSBarry Smith       dk[26]+= uik[12]*u24+ uik[13]*u25+ uik[14]*u26+ uik[15]*u27+ uik[16]*u28+ uik[17]*u29;
1791b3064deSBarry Smith       dk[27]+= uik[18]*u24+ uik[19]*u25+ uik[20]*u26+ uik[21]*u27+ uik[22]*u28+ uik[23]*u29;
1801b3064deSBarry Smith       dk[28]+= uik[24]*u24+ uik[25]*u25+ uik[26]*u26+ uik[27]*u27+ uik[28]*u28+ uik[29]*u29;
1811b3064deSBarry Smith       dk[29]+= uik[30]*u24+ uik[31]*u25+ uik[32]*u26+ uik[33]*u27+ uik[34]*u28+ uik[35]*u29;
18281278733SSatish Balay 
1831b3064deSBarry Smith       dk[30]+=  uik[0]*u30 + uik[1]*u31 + uik[2]*u32 + uik[3]*u33 + uik[4]*u34 + uik[5]*u35;
1841b3064deSBarry Smith       dk[31]+=  uik[6]*u30 + uik[7]*u31 + uik[8]*u32 + uik[9]*u33+ uik[10]*u34+ uik[11]*u35;
1851b3064deSBarry Smith       dk[32]+= uik[12]*u30+ uik[13]*u31+ uik[14]*u32+ uik[15]*u33+ uik[16]*u34+ uik[17]*u35;
1861b3064deSBarry Smith       dk[33]+= uik[18]*u30+ uik[19]*u31+ uik[20]*u32+ uik[21]*u33+ uik[22]*u34+ uik[23]*u35;
1871b3064deSBarry Smith       dk[34]+= uik[24]*u30+ uik[25]*u31+ uik[26]*u32+ uik[27]*u33+ uik[28]*u34+ uik[29]*u35;
1881b3064deSBarry Smith       dk[35]+= uik[30]*u30+ uik[31]*u31+ uik[32]*u32+ uik[33]*u33+ uik[34]*u34+ uik[35]*u35;
18981278733SSatish Balay 
190187a9f4bSHong Zhang       ierr = PetscLogFlops(216*4);CHKERRQ(ierr);
191187a9f4bSHong Zhang 
19281278733SSatish Balay       /* update -U(i,k) */
19381278733SSatish Balay       ierr = PetscMemcpy(ba+ili*36,uik,36*sizeof(MatScalar));CHKERRQ(ierr);
19481278733SSatish Balay 
19581278733SSatish Balay       /* add multiple of row i to k-th row ... */
19681278733SSatish Balay       jmin = ili + 1; jmax = bi[i+1];
19781278733SSatish Balay       if (jmin < jmax){
19881278733SSatish Balay         for (j=jmin; j<jmax; j++) {
19981278733SSatish Balay           /* w += -U(i,k)^T * U_bar(i,j) */
20081278733SSatish Balay           wp = w + bj[j]*36;
20181278733SSatish Balay           u = ba + j*36;
20281278733SSatish Balay 
2031b3064deSBarry Smith 	  u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6];
2041b3064deSBarry Smith 	  u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13];
2051b3064deSBarry Smith 	  u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20];
2061b3064deSBarry Smith 	  u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27];
2071b3064deSBarry Smith 	  u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34];
2081b3064deSBarry Smith 	  u35 = u[35];
20981278733SSatish Balay 
2101b3064deSBarry Smith           wp[0] +=  uik[0]*u0 + uik[1]*u1 + uik[2]*u2 + uik[3]*u3 + uik[4]*u4 + uik[5]*u5;
2111b3064deSBarry Smith           wp[1] +=  uik[6]*u0 + uik[7]*u1 + uik[8]*u2 + uik[9]*u3+ uik[10]*u4+ uik[11]*u5;
2121b3064deSBarry Smith           wp[2] += uik[12]*u0+ uik[13]*u1+ uik[14]*u2+ uik[15]*u3+ uik[16]*u4+ uik[17]*u5;
2131b3064deSBarry Smith           wp[3] += uik[18]*u0+ uik[19]*u1+ uik[20]*u2+ uik[21]*u3+ uik[22]*u4+ uik[23]*u5;
2141b3064deSBarry Smith           wp[4] += uik[24]*u0+ uik[25]*u1+ uik[26]*u2+ uik[27]*u3+ uik[28]*u4+ uik[29]*u5;
2151b3064deSBarry Smith           wp[5] += uik[30]*u0+ uik[31]*u1+ uik[32]*u2+ uik[33]*u3+ uik[34]*u4+ uik[35]*u5;
21681278733SSatish Balay 
2171b3064deSBarry Smith           wp[6] +=  uik[0]*u6 + uik[1]*u7 + uik[2]*u8 + uik[3]*u9 + uik[4]*u10 + uik[5]*u11;
2181b3064deSBarry Smith           wp[7] +=  uik[6]*u6 + uik[7]*u7 + uik[8]*u8 + uik[9]*u9+ uik[10]*u10+ uik[11]*u11;
2191b3064deSBarry Smith           wp[8] += uik[12]*u6+ uik[13]*u7+ uik[14]*u8+ uik[15]*u9+ uik[16]*u10+ uik[17]*u11;
2201b3064deSBarry Smith           wp[9] += uik[18]*u6+ uik[19]*u7+ uik[20]*u8+ uik[21]*u9+ uik[22]*u10+ uik[23]*u11;
2211b3064deSBarry Smith           wp[10]+= uik[24]*u6+ uik[25]*u7+ uik[26]*u8+ uik[27]*u9+ uik[28]*u10+ uik[29]*u11;
2221b3064deSBarry Smith           wp[11]+= uik[30]*u6+ uik[31]*u7+ uik[32]*u8+ uik[33]*u9+ uik[34]*u10+ uik[35]*u11;
22381278733SSatish Balay 
2241b3064deSBarry Smith           wp[12]+=  uik[0]*u12 + uik[1]*u13 + uik[2]*u14 + uik[3]*u15 + uik[4]*u16 + uik[5]*u17;
2251b3064deSBarry Smith           wp[13]+=  uik[6]*u12 + uik[7]*u13 + uik[8]*u14 + uik[9]*u15+ uik[10]*u16+ uik[11]*u17;
2261b3064deSBarry Smith           wp[14]+= uik[12]*u12+ uik[13]*u13+ uik[14]*u14+ uik[15]*u15+ uik[16]*u16+ uik[17]*u17;
2271b3064deSBarry Smith           wp[15]+= uik[18]*u12+ uik[19]*u13+ uik[20]*u14+ uik[21]*u15+ uik[22]*u16+ uik[23]*u17;
2281b3064deSBarry Smith           wp[16]+= uik[24]*u12+ uik[25]*u13+ uik[26]*u14+ uik[27]*u15+ uik[28]*u16+ uik[29]*u17;
2291b3064deSBarry Smith           wp[17]+= uik[30]*u12+ uik[31]*u13+ uik[32]*u14+ uik[33]*u15+ uik[34]*u16+ uik[35]*u17;
23081278733SSatish Balay 
2311b3064deSBarry Smith           wp[18]+=  uik[0]*u18 + uik[1]*u19 + uik[2]*u20 + uik[3]*u21 + uik[4]*u22 + uik[5]*u23;
2321b3064deSBarry Smith           wp[19]+=  uik[6]*u18 + uik[7]*u19 + uik[8]*u20 + uik[9]*u21+ uik[10]*u22+ uik[11]*u23;
2331b3064deSBarry Smith           wp[20]+= uik[12]*u18+ uik[13]*u19+ uik[14]*u20+ uik[15]*u21+ uik[16]*u22+ uik[17]*u23;
2341b3064deSBarry Smith           wp[21]+= uik[18]*u18+ uik[19]*u19+ uik[20]*u20+ uik[21]*u21+ uik[22]*u22+ uik[23]*u23;
2351b3064deSBarry Smith           wp[22]+= uik[24]*u18+ uik[25]*u19+ uik[26]*u20+ uik[27]*u21+ uik[28]*u22+ uik[29]*u23;
2361b3064deSBarry Smith           wp[23]+= uik[30]*u18+ uik[31]*u19+ uik[32]*u20+ uik[33]*u21+ uik[34]*u22+ uik[35]*u23;
2371b3064deSBarry Smith 
2381b3064deSBarry Smith           wp[24]+=  uik[0]*u24 + uik[1]*u25 + uik[2]*u26 + uik[3]*u27 + uik[4]*u28 + uik[5]*u29;
2391b3064deSBarry Smith           wp[25]+=  uik[6]*u24 + uik[7]*u25 + uik[8]*u26 + uik[9]*u27+ uik[10]*u28+ uik[11]*u29;
2401b3064deSBarry Smith           wp[26]+= uik[12]*u24+ uik[13]*u25+ uik[14]*u26+ uik[15]*u27+ uik[16]*u28+ uik[17]*u29;
2411b3064deSBarry Smith           wp[27]+= uik[18]*u24+ uik[19]*u25+ uik[20]*u26+ uik[21]*u27+ uik[22]*u28+ uik[23]*u29;
2421b3064deSBarry Smith           wp[28]+= uik[24]*u24+ uik[25]*u25+ uik[26]*u26+ uik[27]*u27+ uik[28]*u28+ uik[29]*u29;
2431b3064deSBarry Smith           wp[29]+= uik[30]*u24+ uik[31]*u25+ uik[32]*u26+ uik[33]*u27+ uik[34]*u28+ uik[35]*u29;
2441b3064deSBarry Smith 
2451b3064deSBarry Smith           wp[30]+=  uik[0]*u30 + uik[1]*u31 + uik[2]*u32 + uik[3]*u33 + uik[4]*u34 + uik[5]*u35;
2461b3064deSBarry Smith           wp[31]+=  uik[6]*u30 + uik[7]*u31 + uik[8]*u32 + uik[9]*u33+ uik[10]*u34+ uik[11]*u35;
2471b3064deSBarry Smith           wp[32]+= uik[12]*u30+ uik[13]*u31+ uik[14]*u32+ uik[15]*u33+ uik[16]*u34+ uik[17]*u35;
2481b3064deSBarry Smith           wp[33]+= uik[18]*u30+ uik[19]*u31+ uik[20]*u32+ uik[21]*u33+ uik[22]*u34+ uik[23]*u35;
2491b3064deSBarry Smith           wp[34]+= uik[24]*u30+ uik[25]*u31+ uik[26]*u32+ uik[27]*u33+ uik[28]*u34+ uik[29]*u35;
2501b3064deSBarry Smith           wp[35]+= uik[30]*u30+ uik[31]*u31+ uik[32]*u32+ uik[33]*u33+ uik[34]*u34+ uik[35]*u35;
25181278733SSatish Balay         }
252187a9f4bSHong Zhang         ierr = PetscLogFlops(2*216*(jmax-jmin));CHKERRQ(ierr);
25381278733SSatish Balay 
25481278733SSatish Balay         /* ... add i to row list for next nonzero entry */
25581278733SSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
25681278733SSatish Balay         j     = bj[jmin];
25781278733SSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
25881278733SSatish Balay       }
25981278733SSatish Balay       i = nexti;
26081278733SSatish Balay     }
26181278733SSatish Balay 
26281278733SSatish Balay     /* save nonzero entries in k-th row of U ... */
26381278733SSatish Balay 
26481278733SSatish Balay     /* invert diagonal block */
26581278733SSatish Balay     d = ba+k*36;
26681278733SSatish Balay     ierr = PetscMemcpy(d,dk,36*sizeof(MatScalar));CHKERRQ(ierr);
26762bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_6(d,shift);CHKERRQ(ierr);
26881278733SSatish Balay 
26981278733SSatish Balay     jmin = bi[k]; jmax = bi[k+1];
27081278733SSatish Balay     if (jmin < jmax) {
27181278733SSatish Balay       for (j=jmin; j<jmax; j++){
27281278733SSatish Balay          vj = bj[j];           /* block col. index of U */
27381278733SSatish Balay          u   = ba + j*36;
27481278733SSatish Balay          wp = w + vj*36;
27581278733SSatish Balay          for (k1=0; k1<36; k1++){
27681278733SSatish Balay            *u++        = *wp;
27781278733SSatish Balay            *wp++ = 0.0;
27881278733SSatish Balay          }
27981278733SSatish Balay       }
28081278733SSatish Balay 
28181278733SSatish Balay       /* ... add k to row list for first nonzero entry in k-th row */
28281278733SSatish Balay       il[k] = jmin;
28381278733SSatish Balay       i     = bj[jmin];
28481278733SSatish Balay       jl[k] = jl[i]; jl[i] = k;
28581278733SSatish Balay     }
28681278733SSatish Balay   }
28781278733SSatish Balay 
28881278733SSatish Balay   ierr = PetscFree(w);CHKERRQ(ierr);
28981278733SSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
29081278733SSatish Balay   ierr = PetscFree(dk);CHKERRQ(ierr);
29181278733SSatish Balay   if (a->permute){
29281278733SSatish Balay     ierr = PetscFree(aa);CHKERRQ(ierr);
29381278733SSatish Balay   }
29481278733SSatish Balay 
29581278733SSatish Balay   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
296db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqSBAIJ_6;
297db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolve_SeqSBAIJ_6;
29881278733SSatish Balay   C->assembled = PETSC_TRUE;
29981278733SSatish Balay   C->preallocated = PETSC_TRUE;
300efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*216*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
30181278733SSatish Balay   PetscFunctionReturn(0);
30281278733SSatish Balay }
303