xref: /petsc/src/mat/impls/aij/seq/inode.c (revision c6b1f4109a937c85e107e5fb3f283dd757d0671d)
14c1414c8SBarry Smith #define PETSCMAT_DLL
24c1414c8SBarry Smith 
34c1414c8SBarry Smith /*
44c1414c8SBarry Smith   This file provides high performance routines for the Inode format (compressed sparse row)
54c1414c8SBarry Smith   by taking advantage of rows with identical nonzero structure (I-nodes).
64c1414c8SBarry Smith */
74c1414c8SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
84c1414c8SBarry Smith 
94c1414c8SBarry Smith #undef __FUNCT__
104c1414c8SBarry Smith #define __FUNCT__ "Mat_CreateColInode"
114c1414c8SBarry Smith static PetscErrorCode Mat_CreateColInode(Mat A,PetscInt* size,PetscInt ** ns)
124c1414c8SBarry Smith {
134c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
144c1414c8SBarry Smith   PetscErrorCode ierr;
154c1414c8SBarry Smith   PetscInt       i,count,m,n,min_mn,*ns_row,*ns_col;
164c1414c8SBarry Smith 
174c1414c8SBarry Smith   PetscFunctionBegin;
184c1414c8SBarry Smith   n      = A->cmap.n;
194c1414c8SBarry Smith   m      = A->rmap.n;
204c1414c8SBarry Smith   ns_row = a->inode.size;
214c1414c8SBarry Smith 
224c1414c8SBarry Smith   min_mn = (m < n) ? m : n;
234c1414c8SBarry Smith   if (!ns) {
244c1414c8SBarry Smith     for (count=0,i=0; count<min_mn; count+=ns_row[i],i++);
254c1414c8SBarry Smith     for(; count+1 < n; count++,i++);
264c1414c8SBarry Smith     if (count < n)  {
274c1414c8SBarry Smith       i++;
284c1414c8SBarry Smith     }
294c1414c8SBarry Smith     *size = i;
304c1414c8SBarry Smith     PetscFunctionReturn(0);
314c1414c8SBarry Smith   }
324c1414c8SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ns_col);CHKERRQ(ierr);
334c1414c8SBarry Smith 
344c1414c8SBarry Smith   /* Use the same row structure wherever feasible. */
354c1414c8SBarry Smith   for (count=0,i=0; count<min_mn; count+=ns_row[i],i++) {
364c1414c8SBarry Smith     ns_col[i] = ns_row[i];
374c1414c8SBarry Smith   }
384c1414c8SBarry Smith 
394c1414c8SBarry Smith   /* if m < n; pad up the remainder with inode_limit */
404c1414c8SBarry Smith   for(; count+1 < n; count++,i++) {
414c1414c8SBarry Smith     ns_col[i] = 1;
424c1414c8SBarry Smith   }
434c1414c8SBarry Smith   /* The last node is the odd ball. padd it up with the remaining rows; */
444c1414c8SBarry Smith   if (count < n)  {
454c1414c8SBarry Smith     ns_col[i] = n - count;
464c1414c8SBarry Smith     i++;
474c1414c8SBarry Smith   } else if (count > n) {
484c1414c8SBarry Smith     /* Adjust for the over estimation */
494c1414c8SBarry Smith     ns_col[i-1] += n - count;
504c1414c8SBarry Smith   }
514c1414c8SBarry Smith   *size = i;
524c1414c8SBarry Smith   *ns   = ns_col;
534c1414c8SBarry Smith   PetscFunctionReturn(0);
544c1414c8SBarry Smith }
554c1414c8SBarry Smith 
564c1414c8SBarry Smith 
574c1414c8SBarry Smith /*
584c1414c8SBarry Smith       This builds symmetric version of nonzero structure,
594c1414c8SBarry Smith */
604c1414c8SBarry Smith #undef __FUNCT__
614c1414c8SBarry Smith #define __FUNCT__ "MatGetRowIJ_Inode_Symmetric"
624c1414c8SBarry Smith static PetscErrorCode MatGetRowIJ_Inode_Symmetric(Mat A,PetscInt *iia[],PetscInt *jja[],PetscInt ishift,PetscInt oshift)
634c1414c8SBarry Smith {
644c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
654c1414c8SBarry Smith   PetscErrorCode ierr;
664c1414c8SBarry Smith   PetscInt       *work,*ia,*ja,*j,nz,nslim_row,nslim_col,m,row,col,*jmax,n;
674c1414c8SBarry Smith   PetscInt       *tns,*tvc,*ns_row = a->inode.size,*ns_col,nsz,i1,i2,*ai= a->i,*aj = a->j;
684c1414c8SBarry Smith 
694c1414c8SBarry Smith   PetscFunctionBegin;
704c1414c8SBarry Smith   nslim_row = a->inode.node_count;
714c1414c8SBarry Smith   m         = A->rmap.n;
724c1414c8SBarry Smith   n         = A->cmap.n;
734c1414c8SBarry Smith   if (m != n) SETERRQ(PETSC_ERR_SUP,"MatGetRowIJ_Inode_Symmetric: Matrix should be square");
744c1414c8SBarry Smith 
754c1414c8SBarry Smith   /* Use the row_inode as column_inode */
764c1414c8SBarry Smith   nslim_col = nslim_row;
774c1414c8SBarry Smith   ns_col    = ns_row;
784c1414c8SBarry Smith 
794c1414c8SBarry Smith   /* allocate space for reformated inode structure */
804c1414c8SBarry Smith   ierr = PetscMalloc((nslim_col+1)*sizeof(PetscInt),&tns);CHKERRQ(ierr);
814c1414c8SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&tvc);CHKERRQ(ierr);
824c1414c8SBarry Smith   for (i1=0,tns[0]=0; i1<nslim_col; ++i1) tns[i1+1] = tns[i1]+ ns_row[i1];
834c1414c8SBarry Smith 
844c1414c8SBarry Smith   for (i1=0,col=0; i1<nslim_col; ++i1){
854c1414c8SBarry Smith     nsz = ns_col[i1];
864c1414c8SBarry Smith     for (i2=0; i2<nsz; ++i2,++col)
874c1414c8SBarry Smith       tvc[col] = i1;
884c1414c8SBarry Smith   }
894c1414c8SBarry Smith   /* allocate space for row pointers */
904c1414c8SBarry Smith   ierr = PetscMalloc((nslim_row+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
914c1414c8SBarry Smith   *iia = ia;
924c1414c8SBarry Smith   ierr = PetscMemzero(ia,(nslim_row+1)*sizeof(PetscInt));CHKERRQ(ierr);
934c1414c8SBarry Smith   ierr = PetscMalloc((nslim_row+1)*sizeof(PetscInt),&work);CHKERRQ(ierr);
944c1414c8SBarry Smith 
954c1414c8SBarry Smith   /* determine the number of columns in each row */
964c1414c8SBarry Smith   ia[0] = oshift;
974c1414c8SBarry Smith   for (i1=0,row=0 ; i1<nslim_row; row+=ns_row[i1],i1++) {
984c1414c8SBarry Smith 
994c1414c8SBarry Smith     j    = aj + ai[row] + ishift;
1004c1414c8SBarry Smith     jmax = aj + ai[row+1] + ishift;
1014c1414c8SBarry Smith     i2   = 0;
1024c1414c8SBarry Smith     col  = *j++ + ishift;
1034c1414c8SBarry Smith     i2   = tvc[col];
1044c1414c8SBarry Smith     while (i2<i1 && j<jmax) { /* 1.[-xx-d-xx--] 2.[-xx-------],off-diagonal elemets */
1054c1414c8SBarry Smith       ia[i1+1]++;
1064c1414c8SBarry Smith       ia[i2+1]++;
1074c1414c8SBarry Smith       i2++;                     /* Start col of next node */
1084c1414c8SBarry Smith       while(((col=*j+ishift)<tns[i2]) && (j<jmax)) ++j;
1094c1414c8SBarry Smith       i2 = tvc[col];
1104c1414c8SBarry Smith     }
1114c1414c8SBarry Smith     if(i2 == i1) ia[i2+1]++;    /* now the diagonal element */
1124c1414c8SBarry Smith   }
1134c1414c8SBarry Smith 
1144c1414c8SBarry Smith   /* shift ia[i] to point to next row */
1154c1414c8SBarry Smith   for (i1=1; i1<nslim_row+1; i1++) {
1164c1414c8SBarry Smith     row        = ia[i1-1];
1174c1414c8SBarry Smith     ia[i1]    += row;
1184c1414c8SBarry Smith     work[i1-1] = row - oshift;
1194c1414c8SBarry Smith   }
1204c1414c8SBarry Smith 
1214c1414c8SBarry Smith   /* allocate space for column pointers */
1224c1414c8SBarry Smith   nz   = ia[nslim_row] + (!ishift);
1234c1414c8SBarry Smith   ierr = PetscMalloc(nz*sizeof(PetscInt),&ja);CHKERRQ(ierr);
1244c1414c8SBarry Smith   *jja = ja;
1254c1414c8SBarry Smith 
1264c1414c8SBarry Smith  /* loop over lower triangular part putting into ja */
1274c1414c8SBarry Smith   for (i1=0,row=0; i1<nslim_row; row += ns_row[i1],i1++) {
1284c1414c8SBarry Smith     j    = aj + ai[row] + ishift;
1294c1414c8SBarry Smith     jmax = aj + ai[row+1] + ishift;
1304c1414c8SBarry Smith     i2   = 0;                     /* Col inode index */
1314c1414c8SBarry Smith     col  = *j++ + ishift;
1324c1414c8SBarry Smith     i2   = tvc[col];
1334c1414c8SBarry Smith     while (i2<i1 && j<jmax) {
1344c1414c8SBarry Smith       ja[work[i2]++] = i1 + oshift;
1354c1414c8SBarry Smith       ja[work[i1]++] = i2 + oshift;
1364c1414c8SBarry Smith       ++i2;
1374c1414c8SBarry Smith       while(((col=*j+ishift)< tns[i2])&&(j<jmax)) ++j; /* Skip rest col indices in this node */
1384c1414c8SBarry Smith       i2 = tvc[col];
1394c1414c8SBarry Smith     }
1404c1414c8SBarry Smith     if (i2 == i1) ja[work[i1]++] = i2 + oshift;
1414c1414c8SBarry Smith 
1424c1414c8SBarry Smith   }
1434c1414c8SBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
1444c1414c8SBarry Smith   ierr = PetscFree(tns);CHKERRQ(ierr);
1454c1414c8SBarry Smith   ierr = PetscFree(tvc);CHKERRQ(ierr);
1464c1414c8SBarry Smith   PetscFunctionReturn(0);
1474c1414c8SBarry Smith }
1484c1414c8SBarry Smith 
1494c1414c8SBarry Smith /*
1504c1414c8SBarry Smith       This builds nonsymmetric version of nonzero structure,
1514c1414c8SBarry Smith */
1524c1414c8SBarry Smith #undef __FUNCT__
1534c1414c8SBarry Smith #define __FUNCT__ "MatGetRowIJ_Inode_Nonsymmetric"
1544c1414c8SBarry Smith static PetscErrorCode MatGetRowIJ_Inode_Nonsymmetric(Mat A,PetscInt *iia[],PetscInt *jja[],PetscInt ishift,PetscInt oshift)
1554c1414c8SBarry Smith {
1564c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1574c1414c8SBarry Smith   PetscErrorCode ierr;
1584c1414c8SBarry Smith   PetscInt       *work,*ia,*ja,*j,nz,nslim_row,n,row,col,*ns_col,nslim_col;
1594c1414c8SBarry Smith   PetscInt       *tns,*tvc,*ns_row = a->inode.size,nsz,i1,i2,*ai= a->i,*aj = a->j;
1604c1414c8SBarry Smith 
1614c1414c8SBarry Smith   PetscFunctionBegin;
1624c1414c8SBarry Smith   nslim_row = a->inode.node_count;
1634c1414c8SBarry Smith   n         = A->cmap.n;
1644c1414c8SBarry Smith 
1654c1414c8SBarry Smith   /* Create The column_inode for this matrix */
1664c1414c8SBarry Smith   ierr = Mat_CreateColInode(A,&nslim_col,&ns_col);CHKERRQ(ierr);
1674c1414c8SBarry Smith 
1684c1414c8SBarry Smith   /* allocate space for reformated column_inode structure */
1694c1414c8SBarry Smith   ierr = PetscMalloc((nslim_col +1)*sizeof(PetscInt),&tns);CHKERRQ(ierr);
1704c1414c8SBarry Smith   ierr = PetscMalloc((n +1)*sizeof(PetscInt),&tvc);CHKERRQ(ierr);
1714c1414c8SBarry Smith   for (i1=0,tns[0]=0; i1<nslim_col; ++i1) tns[i1+1] = tns[i1] + ns_col[i1];
1724c1414c8SBarry Smith 
1734c1414c8SBarry Smith   for (i1=0,col=0; i1<nslim_col; ++i1){
1744c1414c8SBarry Smith     nsz = ns_col[i1];
1754c1414c8SBarry Smith     for (i2=0; i2<nsz; ++i2,++col)
1764c1414c8SBarry Smith       tvc[col] = i1;
1774c1414c8SBarry Smith   }
1784c1414c8SBarry Smith   /* allocate space for row pointers */
1794c1414c8SBarry Smith   ierr = PetscMalloc((nslim_row+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
1804c1414c8SBarry Smith   *iia = ia;
1814c1414c8SBarry Smith   ierr = PetscMemzero(ia,(nslim_row+1)*sizeof(PetscInt));CHKERRQ(ierr);
1824c1414c8SBarry Smith   ierr = PetscMalloc((nslim_row+1)*sizeof(PetscInt),&work);CHKERRQ(ierr);
1834c1414c8SBarry Smith 
1844c1414c8SBarry Smith   /* determine the number of columns in each row */
1854c1414c8SBarry Smith   ia[0] = oshift;
1864c1414c8SBarry Smith   for (i1=0,row=0; i1<nslim_row; row+=ns_row[i1],i1++) {
1874c1414c8SBarry Smith     j   = aj + ai[row] + ishift;
1884c1414c8SBarry Smith     col = *j++ + ishift;
1894c1414c8SBarry Smith     i2  = tvc[col];
1904c1414c8SBarry Smith     nz  = ai[row+1] - ai[row];
1914c1414c8SBarry Smith     while (nz-- > 0) {           /* off-diagonal elemets */
1924c1414c8SBarry Smith       ia[i1+1]++;
1934c1414c8SBarry Smith       i2++;                     /* Start col of next node */
1944c1414c8SBarry Smith       while (((col = *j++ + ishift) < tns[i2]) && nz > 0) {nz--;}
1954c1414c8SBarry Smith       if (nz > 0) i2 = tvc[col];
1964c1414c8SBarry Smith     }
1974c1414c8SBarry Smith   }
1984c1414c8SBarry Smith 
1994c1414c8SBarry Smith   /* shift ia[i] to point to next row */
2004c1414c8SBarry Smith   for (i1=1; i1<nslim_row+1; i1++) {
2014c1414c8SBarry Smith     row        = ia[i1-1];
2024c1414c8SBarry Smith     ia[i1]    += row;
2034c1414c8SBarry Smith     work[i1-1] = row - oshift;
2044c1414c8SBarry Smith   }
2054c1414c8SBarry Smith 
2064c1414c8SBarry Smith   /* allocate space for column pointers */
2074c1414c8SBarry Smith   nz   = ia[nslim_row] + (!ishift);
2084c1414c8SBarry Smith   ierr = PetscMalloc(nz*sizeof(PetscInt),&ja);CHKERRQ(ierr);
2094c1414c8SBarry Smith   *jja = ja;
2104c1414c8SBarry Smith 
2114c1414c8SBarry Smith  /* loop over matrix putting into ja */
2124c1414c8SBarry Smith   for (i1=0,row=0; i1<nslim_row; row+=ns_row[i1],i1++) {
2134c1414c8SBarry Smith     j   = aj + ai[row] + ishift;
2144c1414c8SBarry Smith     i2  = 0;                     /* Col inode index */
2154c1414c8SBarry Smith     col = *j++ + ishift;
2164c1414c8SBarry Smith     i2  = tvc[col];
2174c1414c8SBarry Smith     nz  = ai[row+1] - ai[row];
2184c1414c8SBarry Smith     while (nz-- > 0) {
2194c1414c8SBarry Smith       ja[work[i1]++] = i2 + oshift;
2204c1414c8SBarry Smith       ++i2;
2214c1414c8SBarry Smith       while(((col = *j++ + ishift) < tns[i2]) && nz > 0) {nz--;}
2224c1414c8SBarry Smith       if (nz > 0) i2 = tvc[col];
2234c1414c8SBarry Smith     }
2244c1414c8SBarry Smith   }
2254c1414c8SBarry Smith   ierr = PetscFree(ns_col);CHKERRQ(ierr);
2264c1414c8SBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
2274c1414c8SBarry Smith   ierr = PetscFree(tns);CHKERRQ(ierr);
2284c1414c8SBarry Smith   ierr = PetscFree(tvc);CHKERRQ(ierr);
2294c1414c8SBarry Smith   PetscFunctionReturn(0);
2304c1414c8SBarry Smith }
2314c1414c8SBarry Smith 
2324c1414c8SBarry Smith #undef __FUNCT__
2334c1414c8SBarry Smith #define __FUNCT__ "MatGetRowIJ_Inode"
2344c1414c8SBarry Smith static PetscErrorCode MatGetRowIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
2354c1414c8SBarry Smith {
2364c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
2374c1414c8SBarry Smith   PetscErrorCode ierr;
2384c1414c8SBarry Smith 
2394c1414c8SBarry Smith   PetscFunctionBegin;
2404c1414c8SBarry Smith   *n     = a->inode.node_count;
2414c1414c8SBarry Smith   if (!ia) PetscFunctionReturn(0);
2424c1414c8SBarry Smith 
2434c1414c8SBarry Smith   if (symmetric) {
2444c1414c8SBarry Smith     ierr = MatGetRowIJ_Inode_Symmetric(A,ia,ja,0,oshift);CHKERRQ(ierr);
2454c1414c8SBarry Smith   } else {
2464c1414c8SBarry Smith     ierr = MatGetRowIJ_Inode_Nonsymmetric(A,ia,ja,0,oshift);CHKERRQ(ierr);
2474c1414c8SBarry Smith   }
2484c1414c8SBarry Smith   PetscFunctionReturn(0);
2494c1414c8SBarry Smith }
2504c1414c8SBarry Smith 
2514c1414c8SBarry Smith #undef __FUNCT__
2524c1414c8SBarry Smith #define __FUNCT__ "MatRestoreRowIJ_Inode"
2534c1414c8SBarry Smith static PetscErrorCode MatRestoreRowIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
2544c1414c8SBarry Smith {
2554c1414c8SBarry Smith   PetscErrorCode ierr;
2564c1414c8SBarry Smith 
2574c1414c8SBarry Smith   PetscFunctionBegin;
2584c1414c8SBarry Smith   if (!ia) PetscFunctionReturn(0);
2594c1414c8SBarry Smith   ierr = PetscFree(*ia);CHKERRQ(ierr);
2604c1414c8SBarry Smith   ierr = PetscFree(*ja);CHKERRQ(ierr);
2614c1414c8SBarry Smith   PetscFunctionReturn(0);
2624c1414c8SBarry Smith }
2634c1414c8SBarry Smith 
2644c1414c8SBarry Smith /* ----------------------------------------------------------- */
2654c1414c8SBarry Smith 
2664c1414c8SBarry Smith #undef __FUNCT__
2674c1414c8SBarry Smith #define __FUNCT__ "MatGetColumnIJ_Inode_Nonsymmetric"
2684c1414c8SBarry Smith static PetscErrorCode MatGetColumnIJ_Inode_Nonsymmetric(Mat A,PetscInt *iia[],PetscInt *jja[],PetscInt ishift,PetscInt oshift)
2694c1414c8SBarry Smith {
2704c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
2714c1414c8SBarry Smith   PetscErrorCode ierr;
2724c1414c8SBarry Smith   PetscInt       *work,*ia,*ja,*j,nz,nslim_row, n,row,col,*ns_col,nslim_col;
2734c1414c8SBarry Smith   PetscInt       *tns,*tvc,*ns_row = a->inode.size,nsz,i1,i2,*ai= a->i,*aj = a->j;
2744c1414c8SBarry Smith 
2754c1414c8SBarry Smith   PetscFunctionBegin;
2764c1414c8SBarry Smith   nslim_row = a->inode.node_count;
2774c1414c8SBarry Smith   n         = A->cmap.n;
2784c1414c8SBarry Smith 
2794c1414c8SBarry Smith   /* Create The column_inode for this matrix */
2804c1414c8SBarry Smith   ierr = Mat_CreateColInode(A,&nslim_col,&ns_col);CHKERRQ(ierr);
2814c1414c8SBarry Smith 
2824c1414c8SBarry Smith   /* allocate space for reformated column_inode structure */
2834c1414c8SBarry Smith   ierr = PetscMalloc((nslim_col + 1)*sizeof(PetscInt),&tns);CHKERRQ(ierr);
2844c1414c8SBarry Smith   ierr = PetscMalloc((n + 1)*sizeof(PetscInt),&tvc);CHKERRQ(ierr);
2854c1414c8SBarry Smith   for (i1=0,tns[0]=0; i1<nslim_col; ++i1) tns[i1+1] = tns[i1] + ns_col[i1];
2864c1414c8SBarry Smith 
2874c1414c8SBarry Smith   for (i1=0,col=0; i1<nslim_col; ++i1){
2884c1414c8SBarry Smith     nsz = ns_col[i1];
2894c1414c8SBarry Smith     for (i2=0; i2<nsz; ++i2,++col)
2904c1414c8SBarry Smith       tvc[col] = i1;
2914c1414c8SBarry Smith   }
2924c1414c8SBarry Smith   /* allocate space for column pointers */
2934c1414c8SBarry Smith   ierr = PetscMalloc((nslim_col+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
2944c1414c8SBarry Smith   *iia = ia;
2954c1414c8SBarry Smith   ierr = PetscMemzero(ia,(nslim_col+1)*sizeof(PetscInt));CHKERRQ(ierr);
2964c1414c8SBarry Smith   ierr = PetscMalloc((nslim_col+1)*sizeof(PetscInt),&work);CHKERRQ(ierr);
2974c1414c8SBarry Smith 
2984c1414c8SBarry Smith   /* determine the number of columns in each row */
2994c1414c8SBarry Smith   ia[0] = oshift;
3004c1414c8SBarry Smith   for (i1=0,row=0; i1<nslim_row; row+=ns_row[i1],i1++) {
3014c1414c8SBarry Smith     j   = aj + ai[row] + ishift;
3024c1414c8SBarry Smith     col = *j++ + ishift;
3034c1414c8SBarry Smith     i2  = tvc[col];
3044c1414c8SBarry Smith     nz  = ai[row+1] - ai[row];
3054c1414c8SBarry Smith     while (nz-- > 0) {           /* off-diagonal elemets */
3064c1414c8SBarry Smith       /* ia[i1+1]++; */
3074c1414c8SBarry Smith       ia[i2+1]++;
3084c1414c8SBarry Smith       i2++;
3094c1414c8SBarry Smith       while (((col = *j++ + ishift) < tns[i2]) && nz > 0) {nz--;}
3104c1414c8SBarry Smith       if (nz > 0) i2 = tvc[col];
3114c1414c8SBarry Smith     }
3124c1414c8SBarry Smith   }
3134c1414c8SBarry Smith 
3144c1414c8SBarry Smith   /* shift ia[i] to point to next col */
3154c1414c8SBarry Smith   for (i1=1; i1<nslim_col+1; i1++) {
3164c1414c8SBarry Smith     col        = ia[i1-1];
3174c1414c8SBarry Smith     ia[i1]    += col;
3184c1414c8SBarry Smith     work[i1-1] = col - oshift;
3194c1414c8SBarry Smith   }
3204c1414c8SBarry Smith 
3214c1414c8SBarry Smith   /* allocate space for column pointers */
3224c1414c8SBarry Smith   nz   = ia[nslim_col] + (!ishift);
3234c1414c8SBarry Smith   ierr = PetscMalloc(nz*sizeof(PetscInt),&ja);CHKERRQ(ierr);
3244c1414c8SBarry Smith   *jja = ja;
3254c1414c8SBarry Smith 
3264c1414c8SBarry Smith  /* loop over matrix putting into ja */
3274c1414c8SBarry Smith   for (i1=0,row=0; i1<nslim_row; row+=ns_row[i1],i1++) {
3284c1414c8SBarry Smith     j   = aj + ai[row] + ishift;
3294c1414c8SBarry Smith     i2  = 0;                     /* Col inode index */
3304c1414c8SBarry Smith     col = *j++ + ishift;
3314c1414c8SBarry Smith     i2  = tvc[col];
3324c1414c8SBarry Smith     nz  = ai[row+1] - ai[row];
3334c1414c8SBarry Smith     while (nz-- > 0) {
3344c1414c8SBarry Smith       /* ja[work[i1]++] = i2 + oshift; */
3354c1414c8SBarry Smith       ja[work[i2]++] = i1 + oshift;
3364c1414c8SBarry Smith       i2++;
3374c1414c8SBarry Smith       while(((col = *j++ + ishift) < tns[i2]) && nz > 0) {nz--;}
3384c1414c8SBarry Smith       if (nz > 0) i2 = tvc[col];
3394c1414c8SBarry Smith     }
3404c1414c8SBarry Smith   }
3414c1414c8SBarry Smith   ierr = PetscFree(ns_col);CHKERRQ(ierr);
3424c1414c8SBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
3434c1414c8SBarry Smith   ierr = PetscFree(tns);CHKERRQ(ierr);
3444c1414c8SBarry Smith   ierr = PetscFree(tvc);CHKERRQ(ierr);
3454c1414c8SBarry Smith   PetscFunctionReturn(0);
3464c1414c8SBarry Smith }
3474c1414c8SBarry Smith 
3484c1414c8SBarry Smith #undef __FUNCT__
3494c1414c8SBarry Smith #define __FUNCT__ "MatGetColumnIJ_Inode"
3504c1414c8SBarry Smith static PetscErrorCode MatGetColumnIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
3514c1414c8SBarry Smith {
3524c1414c8SBarry Smith   PetscErrorCode ierr;
3534c1414c8SBarry Smith 
3544c1414c8SBarry Smith   PetscFunctionBegin;
3554c1414c8SBarry Smith   ierr = Mat_CreateColInode(A,n,PETSC_NULL);CHKERRQ(ierr);
3564c1414c8SBarry Smith   if (!ia) PetscFunctionReturn(0);
3574c1414c8SBarry Smith 
3584c1414c8SBarry Smith   if (symmetric) {
3594c1414c8SBarry Smith     /* Since the indices are symmetric it does'nt matter */
3604c1414c8SBarry Smith     ierr = MatGetRowIJ_Inode_Symmetric(A,ia,ja,0,oshift);CHKERRQ(ierr);
3614c1414c8SBarry Smith   } else {
3624c1414c8SBarry Smith     ierr = MatGetColumnIJ_Inode_Nonsymmetric(A,ia,ja,0,oshift);CHKERRQ(ierr);
3634c1414c8SBarry Smith   }
3644c1414c8SBarry Smith   PetscFunctionReturn(0);
3654c1414c8SBarry Smith }
3664c1414c8SBarry Smith 
3674c1414c8SBarry Smith #undef __FUNCT__
3684c1414c8SBarry Smith #define __FUNCT__ "MatRestoreColumnIJ_Inode"
3694c1414c8SBarry Smith static PetscErrorCode MatRestoreColumnIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
3704c1414c8SBarry Smith {
3714c1414c8SBarry Smith   PetscErrorCode ierr;
3724c1414c8SBarry Smith 
3734c1414c8SBarry Smith   PetscFunctionBegin;
3744c1414c8SBarry Smith   if (!ia) PetscFunctionReturn(0);
3754c1414c8SBarry Smith   ierr = PetscFree(*ia);CHKERRQ(ierr);
3764c1414c8SBarry Smith   ierr = PetscFree(*ja);CHKERRQ(ierr);
3774c1414c8SBarry Smith   PetscFunctionReturn(0);
3784c1414c8SBarry Smith }
3794c1414c8SBarry Smith 
3804c1414c8SBarry Smith /* ----------------------------------------------------------- */
3814c1414c8SBarry Smith 
3824c1414c8SBarry Smith #undef __FUNCT__
3834c1414c8SBarry Smith #define __FUNCT__ "MatMult_Inode"
3844c1414c8SBarry Smith static PetscErrorCode MatMult_Inode(Mat A,Vec xx,Vec yy)
3854c1414c8SBarry Smith {
3864c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3874c1414c8SBarry Smith   PetscScalar    sum1,sum2,sum3,sum4,sum5,tmp0,tmp1;
3884c1414c8SBarry Smith   PetscScalar    *v1,*v2,*v3,*v4,*v5,*x,*y;
3894c1414c8SBarry Smith   PetscErrorCode ierr;
3904c1414c8SBarry Smith   PetscInt       *idx,i1,i2,n,i,row,node_max,*ns,*ii,nsz,sz;
3914c1414c8SBarry Smith 
3924c1414c8SBarry Smith #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
3934c1414c8SBarry Smith #pragma disjoint(*x,*y,*v1,*v2,*v3,*v4,*v5)
3944c1414c8SBarry Smith #endif
3954c1414c8SBarry Smith 
3964c1414c8SBarry Smith   PetscFunctionBegin;
3974c1414c8SBarry Smith   if (!a->inode.size) SETERRQ(PETSC_ERR_COR,"Missing Inode Structure");
3984c1414c8SBarry Smith   node_max = a->inode.node_count;
3994c1414c8SBarry Smith   ns       = a->inode.size;     /* Node Size array */
4004c1414c8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4014c1414c8SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
4024c1414c8SBarry Smith   idx  = a->j;
4034c1414c8SBarry Smith   v1   = a->a;
4044c1414c8SBarry Smith   ii   = a->i;
4054c1414c8SBarry Smith 
4064c1414c8SBarry Smith   for (i = 0,row = 0; i< node_max; ++i){
4074c1414c8SBarry Smith     nsz  = ns[i];
4084c1414c8SBarry Smith     n    = ii[1] - ii[0];
4094c1414c8SBarry Smith     ii  += nsz;
4104c1414c8SBarry Smith     sz   = n;                   /* No of non zeros in this row */
4114c1414c8SBarry Smith                                 /* Switch on the size of Node */
4124c1414c8SBarry Smith     switch (nsz){               /* Each loop in 'case' is unrolled */
4134c1414c8SBarry Smith     case 1 :
4144c1414c8SBarry Smith       sum1  = 0;
4154c1414c8SBarry Smith 
4164c1414c8SBarry Smith       for(n = 0; n< sz-1; n+=2) {
4174c1414c8SBarry Smith         i1   = idx[0];          /* The instructions are ordered to */
4184c1414c8SBarry Smith         i2   = idx[1];          /* make the compiler's job easy */
4194c1414c8SBarry Smith         idx += 2;
4204c1414c8SBarry Smith         tmp0 = x[i1];
4214c1414c8SBarry Smith         tmp1 = x[i2];
4224c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
4234c1414c8SBarry Smith        }
4244c1414c8SBarry Smith 
4254c1414c8SBarry Smith       if (n == sz-1){          /* Take care of the last nonzero  */
4264c1414c8SBarry Smith         tmp0  = x[*idx++];
4274c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
4284c1414c8SBarry Smith       }
4294c1414c8SBarry Smith       y[row++]=sum1;
4304c1414c8SBarry Smith       break;
4314c1414c8SBarry Smith     case 2:
4324c1414c8SBarry Smith       sum1  = 0;
4334c1414c8SBarry Smith       sum2  = 0;
4344c1414c8SBarry Smith       v2    = v1 + n;
4354c1414c8SBarry Smith 
4364c1414c8SBarry Smith       for (n = 0; n< sz-1; n+=2) {
4374c1414c8SBarry Smith         i1   = idx[0];
4384c1414c8SBarry Smith         i2   = idx[1];
4394c1414c8SBarry Smith         idx += 2;
4404c1414c8SBarry Smith         tmp0 = x[i1];
4414c1414c8SBarry Smith         tmp1 = x[i2];
4424c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
4434c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
4444c1414c8SBarry Smith       }
4454c1414c8SBarry Smith       if (n == sz-1){
4464c1414c8SBarry Smith         tmp0  = x[*idx++];
4474c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
4484c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
4494c1414c8SBarry Smith       }
4504c1414c8SBarry Smith       y[row++]=sum1;
4514c1414c8SBarry Smith       y[row++]=sum2;
4524c1414c8SBarry Smith       v1      =v2;              /* Since the next block to be processed starts there*/
4534c1414c8SBarry Smith       idx    +=sz;
4544c1414c8SBarry Smith       break;
4554c1414c8SBarry Smith     case 3:
4564c1414c8SBarry Smith       sum1  = 0;
4574c1414c8SBarry Smith       sum2  = 0;
4584c1414c8SBarry Smith       sum3  = 0;
4594c1414c8SBarry Smith       v2    = v1 + n;
4604c1414c8SBarry Smith       v3    = v2 + n;
4614c1414c8SBarry Smith 
4624c1414c8SBarry Smith       for (n = 0; n< sz-1; n+=2) {
4634c1414c8SBarry Smith         i1   = idx[0];
4644c1414c8SBarry Smith         i2   = idx[1];
4654c1414c8SBarry Smith         idx += 2;
4664c1414c8SBarry Smith         tmp0 = x[i1];
4674c1414c8SBarry Smith         tmp1 = x[i2];
4684c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
4694c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
4704c1414c8SBarry Smith         sum3 += v3[0] * tmp0 + v3[1] * tmp1; v3 += 2;
4714c1414c8SBarry Smith       }
4724c1414c8SBarry Smith       if (n == sz-1){
4734c1414c8SBarry Smith         tmp0  = x[*idx++];
4744c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
4754c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
4764c1414c8SBarry Smith         sum3 += *v3++ * tmp0;
4774c1414c8SBarry Smith       }
4784c1414c8SBarry Smith       y[row++]=sum1;
4794c1414c8SBarry Smith       y[row++]=sum2;
4804c1414c8SBarry Smith       y[row++]=sum3;
4814c1414c8SBarry Smith       v1       =v3;             /* Since the next block to be processed starts there*/
4824c1414c8SBarry Smith       idx     +=2*sz;
4834c1414c8SBarry Smith       break;
4844c1414c8SBarry Smith     case 4:
4854c1414c8SBarry Smith       sum1  = 0;
4864c1414c8SBarry Smith       sum2  = 0;
4874c1414c8SBarry Smith       sum3  = 0;
4884c1414c8SBarry Smith       sum4  = 0;
4894c1414c8SBarry Smith       v2    = v1 + n;
4904c1414c8SBarry Smith       v3    = v2 + n;
4914c1414c8SBarry Smith       v4    = v3 + n;
4924c1414c8SBarry Smith 
4934c1414c8SBarry Smith       for (n = 0; n< sz-1; n+=2) {
4944c1414c8SBarry Smith         i1   = idx[0];
4954c1414c8SBarry Smith         i2   = idx[1];
4964c1414c8SBarry Smith         idx += 2;
4974c1414c8SBarry Smith         tmp0 = x[i1];
4984c1414c8SBarry Smith         tmp1 = x[i2];
4994c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2;
5004c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2;
5014c1414c8SBarry Smith         sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2;
5024c1414c8SBarry Smith         sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2;
5034c1414c8SBarry Smith       }
5044c1414c8SBarry Smith       if (n == sz-1){
5054c1414c8SBarry Smith         tmp0  = x[*idx++];
5064c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
5074c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
5084c1414c8SBarry Smith         sum3 += *v3++ * tmp0;
5094c1414c8SBarry Smith         sum4 += *v4++ * tmp0;
5104c1414c8SBarry Smith       }
5114c1414c8SBarry Smith       y[row++]=sum1;
5124c1414c8SBarry Smith       y[row++]=sum2;
5134c1414c8SBarry Smith       y[row++]=sum3;
5144c1414c8SBarry Smith       y[row++]=sum4;
5154c1414c8SBarry Smith       v1      =v4;              /* Since the next block to be processed starts there*/
5164c1414c8SBarry Smith       idx    +=3*sz;
5174c1414c8SBarry Smith       break;
5184c1414c8SBarry Smith     case 5:
5194c1414c8SBarry Smith       sum1  = 0;
5204c1414c8SBarry Smith       sum2  = 0;
5214c1414c8SBarry Smith       sum3  = 0;
5224c1414c8SBarry Smith       sum4  = 0;
5234c1414c8SBarry Smith       sum5  = 0;
5244c1414c8SBarry Smith       v2    = v1 + n;
5254c1414c8SBarry Smith       v3    = v2 + n;
5264c1414c8SBarry Smith       v4    = v3 + n;
5274c1414c8SBarry Smith       v5    = v4 + n;
5284c1414c8SBarry Smith 
5294c1414c8SBarry Smith       for (n = 0; n<sz-1; n+=2) {
5304c1414c8SBarry Smith         i1   = idx[0];
5314c1414c8SBarry Smith         i2   = idx[1];
5324c1414c8SBarry Smith         idx += 2;
5334c1414c8SBarry Smith         tmp0 = x[i1];
5344c1414c8SBarry Smith         tmp1 = x[i2];
5354c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2;
5364c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2;
5374c1414c8SBarry Smith         sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2;
5384c1414c8SBarry Smith         sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2;
5394c1414c8SBarry Smith         sum5 += v5[0] * tmp0 + v5[1] *tmp1; v5 += 2;
5404c1414c8SBarry Smith       }
5414c1414c8SBarry Smith       if (n == sz-1){
5424c1414c8SBarry Smith         tmp0  = x[*idx++];
5434c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
5444c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
5454c1414c8SBarry Smith         sum3 += *v3++ * tmp0;
5464c1414c8SBarry Smith         sum4 += *v4++ * tmp0;
5474c1414c8SBarry Smith         sum5 += *v5++ * tmp0;
5484c1414c8SBarry Smith       }
5494c1414c8SBarry Smith       y[row++]=sum1;
5504c1414c8SBarry Smith       y[row++]=sum2;
5514c1414c8SBarry Smith       y[row++]=sum3;
5524c1414c8SBarry Smith       y[row++]=sum4;
5534c1414c8SBarry Smith       y[row++]=sum5;
5544c1414c8SBarry Smith       v1      =v5;       /* Since the next block to be processed starts there */
5554c1414c8SBarry Smith       idx    +=4*sz;
5564c1414c8SBarry Smith       break;
5574c1414c8SBarry Smith     default :
5584c1414c8SBarry Smith       SETERRQ(PETSC_ERR_COR,"Node size not yet supported");
5594c1414c8SBarry Smith     }
5604c1414c8SBarry Smith   }
5614c1414c8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5624c1414c8SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
5634c1414c8SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->rmap.n);CHKERRQ(ierr);
5644c1414c8SBarry Smith   PetscFunctionReturn(0);
5654c1414c8SBarry Smith }
5664c1414c8SBarry Smith /* ----------------------------------------------------------- */
5674c1414c8SBarry Smith /* Almost same code as the MatMult_Inode() */
5684c1414c8SBarry Smith #undef __FUNCT__
5694c1414c8SBarry Smith #define __FUNCT__ "MatMultAdd_Inode"
5704c1414c8SBarry Smith static PetscErrorCode MatMultAdd_Inode(Mat A,Vec xx,Vec zz,Vec yy)
5714c1414c8SBarry Smith {
5724c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
5734c1414c8SBarry Smith   PetscScalar    sum1,sum2,sum3,sum4,sum5,tmp0,tmp1;
5744c1414c8SBarry Smith   PetscScalar    *v1,*v2,*v3,*v4,*v5,*x,*y,*z,*zt;
5754c1414c8SBarry Smith   PetscErrorCode ierr;
5764c1414c8SBarry Smith   PetscInt       *idx,i1,i2,n,i,row,node_max,*ns,*ii,nsz,sz;
5774c1414c8SBarry Smith 
5784c1414c8SBarry Smith   PetscFunctionBegin;
5794c1414c8SBarry Smith   if (!a->inode.size) SETERRQ(PETSC_ERR_COR,"Missing Inode Structure");
5804c1414c8SBarry Smith   node_max = a->inode.node_count;
5814c1414c8SBarry Smith   ns       = a->inode.size;     /* Node Size array */
5824c1414c8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5834c1414c8SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5844c1414c8SBarry Smith   if (zz != yy) {
5854c1414c8SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5864c1414c8SBarry Smith   } else {
5874c1414c8SBarry Smith     z = y;
5884c1414c8SBarry Smith   }
5894c1414c8SBarry Smith   zt = z;
5904c1414c8SBarry Smith 
5914c1414c8SBarry Smith   idx  = a->j;
5924c1414c8SBarry Smith   v1   = a->a;
5934c1414c8SBarry Smith   ii   = a->i;
5944c1414c8SBarry Smith 
5954c1414c8SBarry Smith   for (i = 0,row = 0; i< node_max; ++i){
5964c1414c8SBarry Smith     nsz  = ns[i];
5974c1414c8SBarry Smith     n    = ii[1] - ii[0];
5984c1414c8SBarry Smith     ii  += nsz;
5994c1414c8SBarry Smith     sz   = n;                   /* No of non zeros in this row */
6004c1414c8SBarry Smith                                 /* Switch on the size of Node */
6014c1414c8SBarry Smith     switch (nsz){               /* Each loop in 'case' is unrolled */
6024c1414c8SBarry Smith     case 1 :
6034c1414c8SBarry Smith       sum1  = *zt++;
6044c1414c8SBarry Smith 
6054c1414c8SBarry Smith       for(n = 0; n< sz-1; n+=2) {
6064c1414c8SBarry Smith         i1   = idx[0];          /* The instructions are ordered to */
6074c1414c8SBarry Smith         i2   = idx[1];          /* make the compiler's job easy */
6084c1414c8SBarry Smith         idx += 2;
6094c1414c8SBarry Smith         tmp0 = x[i1];
6104c1414c8SBarry Smith         tmp1 = x[i2];
6114c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
6124c1414c8SBarry Smith        }
6134c1414c8SBarry Smith 
6144c1414c8SBarry Smith       if(n   == sz-1){          /* Take care of the last nonzero  */
6154c1414c8SBarry Smith         tmp0  = x[*idx++];
6164c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
6174c1414c8SBarry Smith       }
6184c1414c8SBarry Smith       y[row++]=sum1;
6194c1414c8SBarry Smith       break;
6204c1414c8SBarry Smith     case 2:
6214c1414c8SBarry Smith       sum1  = *zt++;
6224c1414c8SBarry Smith       sum2  = *zt++;
6234c1414c8SBarry Smith       v2    = v1 + n;
6244c1414c8SBarry Smith 
6254c1414c8SBarry Smith       for(n = 0; n< sz-1; n+=2) {
6264c1414c8SBarry Smith         i1   = idx[0];
6274c1414c8SBarry Smith         i2   = idx[1];
6284c1414c8SBarry Smith         idx += 2;
6294c1414c8SBarry Smith         tmp0 = x[i1];
6304c1414c8SBarry Smith         tmp1 = x[i2];
6314c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
6324c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
6334c1414c8SBarry Smith       }
6344c1414c8SBarry Smith       if(n   == sz-1){
6354c1414c8SBarry Smith         tmp0  = x[*idx++];
6364c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
6374c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
6384c1414c8SBarry Smith       }
6394c1414c8SBarry Smith       y[row++]=sum1;
6404c1414c8SBarry Smith       y[row++]=sum2;
6414c1414c8SBarry Smith       v1      =v2;              /* Since the next block to be processed starts there*/
6424c1414c8SBarry Smith       idx    +=sz;
6434c1414c8SBarry Smith       break;
6444c1414c8SBarry Smith     case 3:
6454c1414c8SBarry Smith       sum1  = *zt++;
6464c1414c8SBarry Smith       sum2  = *zt++;
6474c1414c8SBarry Smith       sum3  = *zt++;
6484c1414c8SBarry Smith       v2    = v1 + n;
6494c1414c8SBarry Smith       v3    = v2 + n;
6504c1414c8SBarry Smith 
6514c1414c8SBarry Smith       for (n = 0; n< sz-1; n+=2) {
6524c1414c8SBarry Smith         i1   = idx[0];
6534c1414c8SBarry Smith         i2   = idx[1];
6544c1414c8SBarry Smith         idx += 2;
6554c1414c8SBarry Smith         tmp0 = x[i1];
6564c1414c8SBarry Smith         tmp1 = x[i2];
6574c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
6584c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
6594c1414c8SBarry Smith         sum3 += v3[0] * tmp0 + v3[1] * tmp1; v3 += 2;
6604c1414c8SBarry Smith       }
6614c1414c8SBarry Smith       if (n == sz-1){
6624c1414c8SBarry Smith         tmp0  = x[*idx++];
6634c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
6644c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
6654c1414c8SBarry Smith         sum3 += *v3++ * tmp0;
6664c1414c8SBarry Smith       }
6674c1414c8SBarry Smith       y[row++]=sum1;
6684c1414c8SBarry Smith       y[row++]=sum2;
6694c1414c8SBarry Smith       y[row++]=sum3;
6704c1414c8SBarry Smith       v1       =v3;             /* Since the next block to be processed starts there*/
6714c1414c8SBarry Smith       idx     +=2*sz;
6724c1414c8SBarry Smith       break;
6734c1414c8SBarry Smith     case 4:
6744c1414c8SBarry Smith       sum1  = *zt++;
6754c1414c8SBarry Smith       sum2  = *zt++;
6764c1414c8SBarry Smith       sum3  = *zt++;
6774c1414c8SBarry Smith       sum4  = *zt++;
6784c1414c8SBarry Smith       v2    = v1 + n;
6794c1414c8SBarry Smith       v3    = v2 + n;
6804c1414c8SBarry Smith       v4    = v3 + n;
6814c1414c8SBarry Smith 
6824c1414c8SBarry Smith       for (n = 0; n< sz-1; n+=2) {
6834c1414c8SBarry Smith         i1   = idx[0];
6844c1414c8SBarry Smith         i2   = idx[1];
6854c1414c8SBarry Smith         idx += 2;
6864c1414c8SBarry Smith         tmp0 = x[i1];
6874c1414c8SBarry Smith         tmp1 = x[i2];
6884c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2;
6894c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2;
6904c1414c8SBarry Smith         sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2;
6914c1414c8SBarry Smith         sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2;
6924c1414c8SBarry Smith       }
6934c1414c8SBarry Smith       if (n == sz-1){
6944c1414c8SBarry Smith         tmp0  = x[*idx++];
6954c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
6964c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
6974c1414c8SBarry Smith         sum3 += *v3++ * tmp0;
6984c1414c8SBarry Smith         sum4 += *v4++ * tmp0;
6994c1414c8SBarry Smith       }
7004c1414c8SBarry Smith       y[row++]=sum1;
7014c1414c8SBarry Smith       y[row++]=sum2;
7024c1414c8SBarry Smith       y[row++]=sum3;
7034c1414c8SBarry Smith       y[row++]=sum4;
7044c1414c8SBarry Smith       v1      =v4;              /* Since the next block to be processed starts there*/
7054c1414c8SBarry Smith       idx    +=3*sz;
7064c1414c8SBarry Smith       break;
7074c1414c8SBarry Smith     case 5:
7084c1414c8SBarry Smith       sum1  = *zt++;
7094c1414c8SBarry Smith       sum2  = *zt++;
7104c1414c8SBarry Smith       sum3  = *zt++;
7114c1414c8SBarry Smith       sum4  = *zt++;
7124c1414c8SBarry Smith       sum5  = *zt++;
7134c1414c8SBarry Smith       v2    = v1 + n;
7144c1414c8SBarry Smith       v3    = v2 + n;
7154c1414c8SBarry Smith       v4    = v3 + n;
7164c1414c8SBarry Smith       v5    = v4 + n;
7174c1414c8SBarry Smith 
7184c1414c8SBarry Smith       for (n = 0; n<sz-1; n+=2) {
7194c1414c8SBarry Smith         i1   = idx[0];
7204c1414c8SBarry Smith         i2   = idx[1];
7214c1414c8SBarry Smith         idx += 2;
7224c1414c8SBarry Smith         tmp0 = x[i1];
7234c1414c8SBarry Smith         tmp1 = x[i2];
7244c1414c8SBarry Smith         sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2;
7254c1414c8SBarry Smith         sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2;
7264c1414c8SBarry Smith         sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2;
7274c1414c8SBarry Smith         sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2;
7284c1414c8SBarry Smith         sum5 += v5[0] * tmp0 + v5[1] *tmp1; v5 += 2;
7294c1414c8SBarry Smith       }
7304c1414c8SBarry Smith       if(n   == sz-1){
7314c1414c8SBarry Smith         tmp0  = x[*idx++];
7324c1414c8SBarry Smith         sum1 += *v1++ * tmp0;
7334c1414c8SBarry Smith         sum2 += *v2++ * tmp0;
7344c1414c8SBarry Smith         sum3 += *v3++ * tmp0;
7354c1414c8SBarry Smith         sum4 += *v4++ * tmp0;
7364c1414c8SBarry Smith         sum5 += *v5++ * tmp0;
7374c1414c8SBarry Smith       }
7384c1414c8SBarry Smith       y[row++]=sum1;
7394c1414c8SBarry Smith       y[row++]=sum2;
7404c1414c8SBarry Smith       y[row++]=sum3;
7414c1414c8SBarry Smith       y[row++]=sum4;
7424c1414c8SBarry Smith       y[row++]=sum5;
7434c1414c8SBarry Smith       v1      =v5;       /* Since the next block to be processed starts there */
7444c1414c8SBarry Smith       idx    +=4*sz;
7454c1414c8SBarry Smith       break;
7464c1414c8SBarry Smith     default :
7474c1414c8SBarry Smith       SETERRQ(PETSC_ERR_COR,"Node size not yet supported");
7484c1414c8SBarry Smith     }
7494c1414c8SBarry Smith   }
7504c1414c8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7514c1414c8SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7524c1414c8SBarry Smith   if (zz != yy) {
7534c1414c8SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7544c1414c8SBarry Smith   }
7554c1414c8SBarry Smith   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
7564c1414c8SBarry Smith   PetscFunctionReturn(0);
7574c1414c8SBarry Smith }
7584c1414c8SBarry Smith 
7594c1414c8SBarry Smith /* ----------------------------------------------------------- */
7604c1414c8SBarry Smith #undef __FUNCT__
7614c1414c8SBarry Smith #define __FUNCT__ "MatSolve_Inode"
7624c1414c8SBarry Smith PetscErrorCode MatSolve_Inode(Mat A,Vec bb,Vec xx)
7634c1414c8SBarry Smith {
7644c1414c8SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
7654c1414c8SBarry Smith   IS             iscol = a->col,isrow = a->row;
7664c1414c8SBarry Smith   PetscErrorCode ierr;
7674c1414c8SBarry Smith   PetscInt       *r,*c,i,j,n = A->rmap.n,*ai = a->i,nz,*a_j = a->j;
7684c1414c8SBarry Smith   PetscInt       node_max,*ns,row,nsz,aii,*vi,*ad,*aj,i0,i1,*rout,*cout;
7694c1414c8SBarry Smith   PetscScalar    *x,*b,*a_a = a->a,*tmp,*tmps,*aa,tmp0,tmp1;
7704c1414c8SBarry Smith   PetscScalar    sum1,sum2,sum3,sum4,sum5,*v1,*v2,*v3,*v4,*v5;
7714c1414c8SBarry Smith 
7724c1414c8SBarry Smith   PetscFunctionBegin;
7734c1414c8SBarry Smith   if (A->factor!=FACTOR_LU) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unfactored matrix");
7744c1414c8SBarry Smith   if (!a->inode.size) SETERRQ(PETSC_ERR_COR,"Missing Inode Structure");
7754c1414c8SBarry Smith   node_max = a->inode.node_count;
7764c1414c8SBarry Smith   ns       = a->inode.size;     /* Node Size array */
7774c1414c8SBarry Smith 
7784c1414c8SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7794c1414c8SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
7804c1414c8SBarry Smith   tmp  = a->solve_work;
7814c1414c8SBarry Smith 
7824c1414c8SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7834c1414c8SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
7844c1414c8SBarry Smith 
7854c1414c8SBarry Smith   /* forward solve the lower triangular */
7864c1414c8SBarry Smith   tmps = tmp ;
7874c1414c8SBarry Smith   aa   = a_a ;
7884c1414c8SBarry Smith   aj   = a_j ;
7894c1414c8SBarry Smith   ad   = a->diag;
7904c1414c8SBarry Smith 
7914c1414c8SBarry Smith   for (i = 0,row = 0; i< node_max; ++i){
7924c1414c8SBarry Smith     nsz = ns[i];
7934c1414c8SBarry Smith     aii = ai[row];
7944c1414c8SBarry Smith     v1  = aa + aii;
7954c1414c8SBarry Smith     vi  = aj + aii;
7964c1414c8SBarry Smith     nz  = ad[row]- aii;
7974c1414c8SBarry Smith 
7984c1414c8SBarry Smith     switch (nsz){               /* Each loop in 'case' is unrolled */
7994c1414c8SBarry Smith     case 1 :
8004c1414c8SBarry Smith       sum1 = b[*r++];
8014c1414c8SBarry Smith       /*      while (nz--) sum1 -= *v1++ *tmps[*vi++];*/
8024c1414c8SBarry Smith       for(j=0; j<nz-1; j+=2){
8034c1414c8SBarry Smith         i0   = vi[0];
8044c1414c8SBarry Smith         i1   = vi[1];
8054c1414c8SBarry Smith         vi  +=2;
8064c1414c8SBarry Smith         tmp0 = tmps[i0];
8074c1414c8SBarry Smith         tmp1 = tmps[i1];
8084c1414c8SBarry Smith         sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
8094c1414c8SBarry Smith       }
8104c1414c8SBarry Smith       if(j == nz-1){
8114c1414c8SBarry Smith         tmp0 = tmps[*vi++];
8124c1414c8SBarry Smith         sum1 -= *v1++ *tmp0;
8134c1414c8SBarry Smith       }
8144c1414c8SBarry Smith       tmp[row ++]=sum1;
8154c1414c8SBarry Smith       break;
8164c1414c8SBarry Smith     case 2:
8174c1414c8SBarry Smith       sum1 = b[*r++];
8184c1414c8SBarry Smith       sum2 = b[*r++];
8194c1414c8SBarry Smith       v2   = aa + ai[row+1];
8204c1414c8SBarry Smith 
8214c1414c8SBarry Smith       for(j=0; j<nz-1; j+=2){
8224c1414c8SBarry Smith         i0   = vi[0];
8234c1414c8SBarry Smith         i1   = vi[1];
8244c1414c8SBarry Smith         vi  +=2;
8254c1414c8SBarry Smith         tmp0 = tmps[i0];
8264c1414c8SBarry Smith         tmp1 = tmps[i1];
8274c1414c8SBarry Smith         sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
8284c1414c8SBarry Smith         sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
8294c1414c8SBarry Smith       }
8304c1414c8SBarry Smith       if(j == nz-1){
8314c1414c8SBarry Smith         tmp0 = tmps[*vi++];
8324c1414c8SBarry Smith         sum1 -= *v1++ *tmp0;
8334c1414c8SBarry Smith         sum2 -= *v2++ *tmp0;
8344c1414c8SBarry Smith       }
8354c1414c8SBarry Smith       sum2 -= *v2++ * sum1;
8364c1414c8SBarry Smith       tmp[row ++]=sum1;
8374c1414c8SBarry Smith       tmp[row ++]=sum2;
8384c1414c8SBarry Smith       break;
8394c1414c8SBarry Smith     case 3:
8404c1414c8SBarry Smith       sum1 = b[*r++];
8414c1414c8SBarry Smith       sum2 = b[*r++];
8424c1414c8SBarry Smith       sum3 = b[*r++];
8434c1414c8SBarry Smith       v2   = aa + ai[row+1];
8444c1414c8SBarry Smith       v3   = aa + ai[row+2];
8454c1414c8SBarry Smith 
8464c1414c8SBarry Smith       for (j=0; j<nz-1; j+=2){
8474c1414c8SBarry Smith         i0   = vi[0];
8484c1414c8SBarry Smith         i1   = vi[1];
8494c1414c8SBarry Smith         vi  +=2;
8504c1414c8SBarry Smith         tmp0 = tmps[i0];
8514c1414c8SBarry Smith         tmp1 = tmps[i1];
8524c1414c8SBarry Smith         sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
8534c1414c8SBarry Smith         sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
8544c1414c8SBarry Smith         sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2;
8554c1414c8SBarry Smith       }
8564c1414c8SBarry Smith       if (j == nz-1){
8574c1414c8SBarry Smith         tmp0 = tmps[*vi++];
8584c1414c8SBarry Smith         sum1 -= *v1++ *tmp0;
8594c1414c8SBarry Smith         sum2 -= *v2++ *tmp0;
8604c1414c8SBarry Smith         sum3 -= *v3++ *tmp0;
8614c1414c8SBarry Smith       }
8624c1414c8SBarry Smith       sum2 -= *v2++ * sum1;
8634c1414c8SBarry Smith       sum3 -= *v3++ * sum1;
8644c1414c8SBarry Smith       sum3 -= *v3++ * sum2;
8654c1414c8SBarry Smith       tmp[row ++]=sum1;
8664c1414c8SBarry Smith       tmp[row ++]=sum2;
8674c1414c8SBarry Smith       tmp[row ++]=sum3;
8684c1414c8SBarry Smith       break;
8694c1414c8SBarry Smith 
8704c1414c8SBarry Smith     case 4:
8714c1414c8SBarry Smith       sum1 = b[*r++];
8724c1414c8SBarry Smith       sum2 = b[*r++];
8734c1414c8SBarry Smith       sum3 = b[*r++];
8744c1414c8SBarry Smith       sum4 = b[*r++];
8754c1414c8SBarry Smith       v2   = aa + ai[row+1];
8764c1414c8SBarry Smith       v3   = aa + ai[row+2];
8774c1414c8SBarry Smith       v4   = aa + ai[row+3];
8784c1414c8SBarry Smith 
8794c1414c8SBarry Smith       for (j=0; j<nz-1; j+=2){
8804c1414c8SBarry Smith         i0   = vi[0];
8814c1414c8SBarry Smith         i1   = vi[1];
8824c1414c8SBarry Smith         vi  +=2;
8834c1414c8SBarry Smith         tmp0 = tmps[i0];
8844c1414c8SBarry Smith         tmp1 = tmps[i1];
8854c1414c8SBarry Smith         sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
8864c1414c8SBarry Smith         sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
8874c1414c8SBarry Smith         sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2;
8884c1414c8SBarry Smith         sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2;
8894c1414c8SBarry Smith       }
8904c1414c8SBarry Smith       if (j == nz-1){
8914c1414c8SBarry Smith         tmp0 = tmps[*vi++];
8924c1414c8SBarry Smith         sum1 -= *v1++ *tmp0;
8934c1414c8SBarry Smith         sum2 -= *v2++ *tmp0;
8944c1414c8SBarry Smith         sum3 -= *v3++ *tmp0;
8954c1414c8SBarry Smith         sum4 -= *v4++ *tmp0;
8964c1414c8SBarry Smith       }
8974c1414c8SBarry Smith       sum2 -= *v2++ * sum1;
8984c1414c8SBarry Smith       sum3 -= *v3++ * sum1;
8994c1414c8SBarry Smith       sum4 -= *v4++ * sum1;
9004c1414c8SBarry Smith       sum3 -= *v3++ * sum2;
9014c1414c8SBarry Smith       sum4 -= *v4++ * sum2;
9024c1414c8SBarry Smith       sum4 -= *v4++ * sum3;
9034c1414c8SBarry Smith 
9044c1414c8SBarry Smith       tmp[row ++]=sum1;
9054c1414c8SBarry Smith       tmp[row ++]=sum2;
9064c1414c8SBarry Smith       tmp[row ++]=sum3;
9074c1414c8SBarry Smith       tmp[row ++]=sum4;
9084c1414c8SBarry Smith       break;
9094c1414c8SBarry Smith     case 5:
9104c1414c8SBarry Smith       sum1 = b[*r++];
9114c1414c8SBarry Smith       sum2 = b[*r++];
9124c1414c8SBarry Smith       sum3 = b[*r++];
9134c1414c8SBarry Smith       sum4 = b[*r++];
9144c1414c8SBarry Smith       sum5 = b[*r++];
9154c1414c8SBarry Smith       v2   = aa + ai[row+1];
9164c1414c8SBarry Smith       v3   = aa + ai[row+2];
9174c1414c8SBarry Smith       v4   = aa + ai[row+3];
9184c1414c8SBarry Smith       v5   = aa + ai[row+4];
9194c1414c8SBarry Smith 
9204c1414c8SBarry Smith       for (j=0; j<nz-1; j+=2){
9214c1414c8SBarry Smith         i0   = vi[0];
9224c1414c8SBarry Smith         i1   = vi[1];
9234c1414c8SBarry Smith         vi  +=2;
9244c1414c8SBarry Smith         tmp0 = tmps[i0];
9254c1414c8SBarry Smith         tmp1 = tmps[i1];
9264c1414c8SBarry Smith         sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2;
9274c1414c8SBarry Smith         sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2;
9284c1414c8SBarry Smith         sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2;
9294c1414c8SBarry Smith         sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2;
9304c1414c8SBarry Smith         sum5 -= v5[0] * tmp0 + v5[1] * tmp1; v5 += 2;
9314c1414c8SBarry Smith       }
9324c1414c8SBarry Smith       if (j == nz-1){
9334c1414c8SBarry Smith         tmp0 = tmps[*vi++];
9344c1414c8SBarry Smith         sum1 -= *v1++ *tmp0;
9354c1414c8SBarry Smith         sum2 -= *v2++ *tmp0;
9364c1414c8SBarry Smith         sum3 -= *v3++ *tmp0;
9374c1414c8SBarry Smith         sum4 -= *v4++ *tmp0;
9384c1414c8SBarry Smith         sum5 -= *v5++ *tmp0;
9394c1414c8SBarry Smith       }
9404c1414c8SBarry Smith 
9414c1414c8SBarry Smith       sum2 -= *v2++ * sum1;
9424c1414c8SBarry Smith       sum3 -= *v3++ * sum1;
9434c1414c8SBarry Smith       sum4 -= *v4++ * sum1;
9444c1414c8SBarry Smith       sum5 -= *v5++ * sum1;
9454c1414c8SBarry Smith       sum3 -= *v3++ * sum2;
9464c1414c8SBarry Smith       sum4 -= *v4++ * sum2;
9474c1414c8SBarry Smith       sum5 -= *v5++ * sum2;
9484c1414c8SBarry Smith       sum4 -= *v4++ * sum3;
9494c1414c8SBarry Smith       sum5 -= *v5++ * sum3;
9504c1414c8SBarry Smith       sum5 -= *v5++ * sum4;
9514c1414c8SBarry Smith 
9524c1414c8SBarry Smith       tmp[row ++]=sum1;
9534c1414c8SBarry Smith       tmp[row ++]=sum2;
9544c1414c8SBarry Smith       tmp[row ++]=sum3;
9554c1414c8SBarry Smith       tmp[row ++]=sum4;
9564c1414c8SBarry Smith       tmp[row ++]=sum5;
9574c1414c8SBarry Smith       break;
9584c1414c8SBarry Smith     default:
9594c1414c8SBarry Smith       SETERRQ(PETSC_ERR_COR,"Node size not yet supported \n");
9604c1414c8SBarry Smith     }
9614c1414c8SBarry Smith   }
9624c1414c8SBarry Smith   /* backward solve the upper triangular */
9634c1414c8SBarry Smith   for (i=node_max -1 ,row = n-1 ; i>=0; i--){
9644c1414c8SBarry Smith     nsz = ns[i];
9654c1414c8SBarry Smith     aii = ai[row+1] -1;
9664c1414c8SBarry Smith     v1  = aa + aii;
9674c1414c8SBarry Smith     vi  = aj + aii;
9684c1414c8SBarry Smith     nz  = aii- ad[row];
9694c1414c8SBarry Smith     switch (nsz){               /* Each loop in 'case' is unrolled */
9704c1414c8SBarry Smith     case 1 :
9714c1414c8SBarry Smith       sum1 = tmp[row];
9724c1414c8SBarry Smith 
9734c1414c8SBarry Smith       for(j=nz ; j>1; j-=2){
9744c1414c8SBarry Smith         vi  -=2;
9754c1414c8SBarry Smith         i0   = vi[2];
9764c1414c8SBarry Smith         i1   = vi[1];
9774c1414c8SBarry Smith         tmp0 = tmps[i0];
9784c1414c8SBarry Smith         tmp1 = tmps[i1];
9794c1414c8SBarry Smith         v1   -= 2;
9804c1414c8SBarry Smith         sum1 -= v1[2] * tmp0 + v1[1] * tmp1;
9814c1414c8SBarry Smith       }
9824c1414c8SBarry Smith       if (j==1){
9834c1414c8SBarry Smith         tmp0  = tmps[*vi--];
9844c1414c8SBarry Smith         sum1 -= *v1-- * tmp0;
9854c1414c8SBarry Smith       }
9864c1414c8SBarry Smith       x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--;
9874c1414c8SBarry Smith       break;
9884c1414c8SBarry Smith     case 2 :
9894c1414c8SBarry Smith       sum1 = tmp[row];
9904c1414c8SBarry Smith       sum2 = tmp[row -1];
9914c1414c8SBarry Smith       v2   = aa + ai[row]-1;
9924c1414c8SBarry Smith       for (j=nz ; j>1; j-=2){
9934c1414c8SBarry Smith         vi  -=2;
9944c1414c8SBarry Smith         i0   = vi[2];
9954c1414c8SBarry Smith         i1   = vi[1];
9964c1414c8SBarry Smith         tmp0 = tmps[i0];
9974c1414c8SBarry Smith         tmp1 = tmps[i1];
9984c1414c8SBarry Smith         v1   -= 2;
9994c1414c8SBarry Smith         v2   -= 2;
10004c1414c8SBarry Smith         sum1 -= v1[2] * tmp0 + v1[1] * tmp1;
10014c1414c8SBarry Smith         sum2 -= v2[2] * tmp0 + v2[1] * tmp1;
10024c1414c8SBarry Smith       }
10034c1414c8SBarry Smith       if (j==1){
10044c1414c8SBarry Smith         tmp0  = tmps[*vi--];
10054c1414c8SBarry Smith         sum1 -= *v1-- * tmp0;
10064c1414c8SBarry Smith         sum2 -= *v2-- * tmp0;
10074c1414c8SBarry Smith       }
10084c1414c8SBarry Smith 
10094c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--;
10104c1414c8SBarry Smith       sum2   -= *v2-- * tmp0;
10114c1414c8SBarry Smith       x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--;
10124c1414c8SBarry Smith       break;
10134c1414c8SBarry Smith     case 3 :
10144c1414c8SBarry Smith       sum1 = tmp[row];
10154c1414c8SBarry Smith       sum2 = tmp[row -1];
10164c1414c8SBarry Smith       sum3 = tmp[row -2];
10174c1414c8SBarry Smith       v2   = aa + ai[row]-1;
10184c1414c8SBarry Smith       v3   = aa + ai[row -1]-1;
10194c1414c8SBarry Smith       for (j=nz ; j>1; j-=2){
10204c1414c8SBarry Smith         vi  -=2;
10214c1414c8SBarry Smith         i0   = vi[2];
10224c1414c8SBarry Smith         i1   = vi[1];
10234c1414c8SBarry Smith         tmp0 = tmps[i0];
10244c1414c8SBarry Smith         tmp1 = tmps[i1];
10254c1414c8SBarry Smith         v1   -= 2;
10264c1414c8SBarry Smith         v2   -= 2;
10274c1414c8SBarry Smith         v3   -= 2;
10284c1414c8SBarry Smith         sum1 -= v1[2] * tmp0 + v1[1] * tmp1;
10294c1414c8SBarry Smith         sum2 -= v2[2] * tmp0 + v2[1] * tmp1;
10304c1414c8SBarry Smith         sum3 -= v3[2] * tmp0 + v3[1] * tmp1;
10314c1414c8SBarry Smith       }
10324c1414c8SBarry Smith       if (j==1){
10334c1414c8SBarry Smith         tmp0  = tmps[*vi--];
10344c1414c8SBarry Smith         sum1 -= *v1-- * tmp0;
10354c1414c8SBarry Smith         sum2 -= *v2-- * tmp0;
10364c1414c8SBarry Smith         sum3 -= *v3-- * tmp0;
10374c1414c8SBarry Smith       }
10384c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--;
10394c1414c8SBarry Smith       sum2   -= *v2-- * tmp0;
10404c1414c8SBarry Smith       sum3   -= *v3-- * tmp0;
10414c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--;
10424c1414c8SBarry Smith       sum3   -= *v3-- * tmp0;
10434c1414c8SBarry Smith       x[*c--] = tmp[row] = sum3*a_a[ad[row]]; row--;
10444c1414c8SBarry Smith 
10454c1414c8SBarry Smith       break;
10464c1414c8SBarry Smith     case 4 :
10474c1414c8SBarry Smith       sum1 = tmp[row];
10484c1414c8SBarry Smith       sum2 = tmp[row -1];
10494c1414c8SBarry Smith       sum3 = tmp[row -2];
10504c1414c8SBarry Smith       sum4 = tmp[row -3];
10514c1414c8SBarry Smith       v2   = aa + ai[row]-1;
10524c1414c8SBarry Smith       v3   = aa + ai[row -1]-1;
10534c1414c8SBarry Smith       v4   = aa + ai[row -2]-1;
10544c1414c8SBarry Smith 
10554c1414c8SBarry Smith       for (j=nz ; j>1; j-=2){
10564c1414c8SBarry Smith         vi  -=2;
10574c1414c8SBarry Smith         i0   = vi[2];
10584c1414c8SBarry Smith         i1   = vi[1];
10594c1414c8SBarry Smith         tmp0 = tmps[i0];
10604c1414c8SBarry Smith         tmp1 = tmps[i1];
10614c1414c8SBarry Smith         v1  -= 2;
10624c1414c8SBarry Smith         v2  -= 2;
10634c1414c8SBarry Smith         v3  -= 2;
10644c1414c8SBarry Smith         v4  -= 2;
10654c1414c8SBarry Smith         sum1 -= v1[2] * tmp0 + v1[1] * tmp1;
10664c1414c8SBarry Smith         sum2 -= v2[2] * tmp0 + v2[1] * tmp1;
10674c1414c8SBarry Smith         sum3 -= v3[2] * tmp0 + v3[1] * tmp1;
10684c1414c8SBarry Smith         sum4 -= v4[2] * tmp0 + v4[1] * tmp1;
10694c1414c8SBarry Smith       }
10704c1414c8SBarry Smith       if (j==1){
10714c1414c8SBarry Smith         tmp0  = tmps[*vi--];
10724c1414c8SBarry Smith         sum1 -= *v1-- * tmp0;
10734c1414c8SBarry Smith         sum2 -= *v2-- * tmp0;
10744c1414c8SBarry Smith         sum3 -= *v3-- * tmp0;
10754c1414c8SBarry Smith         sum4 -= *v4-- * tmp0;
10764c1414c8SBarry Smith       }
10774c1414c8SBarry Smith 
10784c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--;
10794c1414c8SBarry Smith       sum2   -= *v2-- * tmp0;
10804c1414c8SBarry Smith       sum3   -= *v3-- * tmp0;
10814c1414c8SBarry Smith       sum4   -= *v4-- * tmp0;
10824c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--;
10834c1414c8SBarry Smith       sum3   -= *v3-- * tmp0;
10844c1414c8SBarry Smith       sum4   -= *v4-- * tmp0;
10854c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum3*a_a[ad[row]]; row--;
10864c1414c8SBarry Smith       sum4   -= *v4-- * tmp0;
10874c1414c8SBarry Smith       x[*c--] = tmp[row] = sum4*a_a[ad[row]]; row--;
10884c1414c8SBarry Smith       break;
10894c1414c8SBarry Smith     case 5 :
10904c1414c8SBarry Smith       sum1 = tmp[row];
10914c1414c8SBarry Smith       sum2 = tmp[row -1];
10924c1414c8SBarry Smith       sum3 = tmp[row -2];
10934c1414c8SBarry Smith       sum4 = tmp[row -3];
10944c1414c8SBarry Smith       sum5 = tmp[row -4];
10954c1414c8SBarry Smith       v2   = aa + ai[row]-1;
10964c1414c8SBarry Smith       v3   = aa + ai[row -1]-1;
10974c1414c8SBarry Smith       v4   = aa + ai[row -2]-1;
10984c1414c8SBarry Smith       v5   = aa + ai[row -3]-1;
10994c1414c8SBarry Smith       for (j=nz ; j>1; j-=2){
11004c1414c8SBarry Smith         vi  -= 2;
11014c1414c8SBarry Smith         i0   = vi[2];
11024c1414c8SBarry Smith         i1   = vi[1];
11034c1414c8SBarry Smith         tmp0 = tmps[i0];
11044c1414c8SBarry Smith         tmp1 = tmps[i1];
11054c1414c8SBarry Smith         v1   -= 2;
11064c1414c8SBarry Smith         v2   -= 2;
11074c1414c8SBarry Smith         v3   -= 2;
11084c1414c8SBarry Smith         v4   -= 2;
11094c1414c8SBarry Smith         v5   -= 2;
11104c1414c8SBarry Smith         sum1 -= v1[2] * tmp0 + v1[1] * tmp1;
11114c1414c8SBarry Smith         sum2 -= v2[2] * tmp0 + v2[1] * tmp1;
11124c1414c8SBarry Smith         sum3 -= v3[2] * tmp0 + v3[1] * tmp1;
11134c1414c8SBarry Smith         sum4 -= v4[2] * tmp0 + v4[1] * tmp1;
11144c1414c8SBarry Smith         sum5 -= v5[2] * tmp0 + v5[1] * tmp1;
11154c1414c8SBarry Smith       }
11164c1414c8SBarry Smith       if (j==1){
11174c1414c8SBarry Smith         tmp0  = tmps[*vi--];
11184c1414c8SBarry Smith         sum1 -= *v1-- * tmp0;
11194c1414c8SBarry Smith         sum2 -= *v2-- * tmp0;
11204c1414c8SBarry Smith         sum3 -= *v3-- * tmp0;
11214c1414c8SBarry Smith         sum4 -= *v4-- * tmp0;
11224c1414c8SBarry Smith         sum5 -= *v5-- * tmp0;
11234c1414c8SBarry Smith       }
11244c1414c8SBarry Smith 
11254c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--;
11264c1414c8SBarry Smith       sum2   -= *v2-- * tmp0;
11274c1414c8SBarry Smith       sum3   -= *v3-- * tmp0;
11284c1414c8SBarry Smith       sum4   -= *v4-- * tmp0;
11294c1414c8SBarry Smith       sum5   -= *v5-- * tmp0;
11304c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--;
11314c1414c8SBarry Smith       sum3   -= *v3-- * tmp0;
11324c1414c8SBarry Smith       sum4   -= *v4-- * tmp0;
11334c1414c8SBarry Smith       sum5   -= *v5-- * tmp0;
11344c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum3*a_a[ad[row]]; row--;
11354c1414c8SBarry Smith       sum4   -= *v4-- * tmp0;
11364c1414c8SBarry Smith       sum5   -= *v5-- * tmp0;
11374c1414c8SBarry Smith       tmp0    = x[*c--] = tmp[row] = sum4*a_a[ad[row]]; row--;
11384c1414c8SBarry Smith       sum5   -= *v5-- * tmp0;
11394c1414c8SBarry Smith       x[*c--] = tmp[row] = sum5*a_a[ad[row]]; row--;
11404c1414c8SBarry Smith       break;
11414c1414c8SBarry Smith     default:
11424c1414c8SBarry Smith       SETERRQ(PETSC_ERR_COR,"Node size not yet supported \n");
11434c1414c8SBarry Smith     }
11444c1414c8SBarry Smith   }
11454c1414c8SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11464c1414c8SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11474c1414c8SBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11484c1414c8SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11494c1414c8SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
11504c1414c8SBarry Smith   PetscFunctionReturn(0);
11514c1414c8SBarry Smith }
11524c1414c8SBarry Smith 
11534c1414c8SBarry Smith #undef __FUNCT__
11544c1414c8SBarry Smith #define __FUNCT__ "MatLUFactorNumeric_Inode"
11554c1414c8SBarry Smith PetscErrorCode MatLUFactorNumeric_Inode(Mat A,MatFactorInfo *info,Mat *B)
11564c1414c8SBarry Smith {
11574c1414c8SBarry Smith   Mat            C = *B;
11584c1414c8SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)C->data;
11594c1414c8SBarry Smith   IS             iscol = b->col,isrow = b->row,isicol = b->icol;
11604c1414c8SBarry Smith   PetscErrorCode ierr;
11614c1414c8SBarry Smith   PetscInt       *r,*ic,*c,n = A->rmap.n,*bi = b->i;
11624c1414c8SBarry Smith   PetscInt       *bj = b->j,*nbj=b->j +1,*ajtmp,*bjtmp,nz,nz_tmp,row,prow;
11634c1414c8SBarry Smith   PetscInt       *ics,i,j,idx,*ai = a->i,*aj = a->j,*bd = b->diag,node_max,nodesz;
11644c1414c8SBarry Smith   PetscInt       *ns,*tmp_vec1,*tmp_vec2,*nsmap,*pj;
11654c1414c8SBarry Smith   PetscScalar    *rtmp1,*rtmp2,*rtmp3,*v1,*v2,*v3,*pc1,*pc2,*pc3,mul1,mul2,mul3;
11664c1414c8SBarry Smith   PetscScalar    tmp,*ba = b->a,*aa = a->a,*pv;
11674c1414c8SBarry Smith   PetscReal      rs=0.0;
11684c1414c8SBarry Smith   LUShift_Ctx    sctx;
11694c1414c8SBarry Smith   PetscInt       newshift;
11704c1414c8SBarry Smith 
11714c1414c8SBarry Smith   PetscFunctionBegin;
11724c1414c8SBarry Smith   sctx.shift_top  = 0;
11734c1414c8SBarry Smith   sctx.nshift_max = 0;
11744c1414c8SBarry Smith   sctx.shift_lo   = 0;
11754c1414c8SBarry Smith   sctx.shift_hi   = 0;
11764c1414c8SBarry Smith 
11774c1414c8SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
11784c1414c8SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
11794c1414c8SBarry Smith   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
11804c1414c8SBarry Smith     sctx.shift_top = 0;
11814c1414c8SBarry Smith     for (i=0; i<n; i++) {
11824c1414c8SBarry Smith       /* calculate rs = sum(|aij|)-RealPart(aii), amt of shift needed for this row */
11834c1414c8SBarry Smith       rs    = 0.0;
11844c1414c8SBarry Smith       ajtmp = aj + ai[i];
11854c1414c8SBarry Smith       rtmp1 = aa + ai[i];
11864c1414c8SBarry Smith       nz = ai[i+1] - ai[i];
11874c1414c8SBarry Smith       for (j=0; j<nz; j++){
11884c1414c8SBarry Smith         if (*ajtmp != i){
11894c1414c8SBarry Smith           rs += PetscAbsScalar(*rtmp1++);
11904c1414c8SBarry Smith         } else {
11914c1414c8SBarry Smith           rs -= PetscRealPart(*rtmp1++);
11924c1414c8SBarry Smith         }
11934c1414c8SBarry Smith         ajtmp++;
11944c1414c8SBarry Smith       }
11954c1414c8SBarry Smith       if (rs>sctx.shift_top) sctx.shift_top = rs;
11964c1414c8SBarry Smith     }
11974c1414c8SBarry Smith     if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12;
11984c1414c8SBarry Smith     sctx.shift_top *= 1.1;
11994c1414c8SBarry Smith     sctx.nshift_max = 5;
12004c1414c8SBarry Smith     sctx.shift_lo   = 0.;
12014c1414c8SBarry Smith     sctx.shift_hi   = 1.;
12024c1414c8SBarry Smith   }
12034c1414c8SBarry Smith   sctx.shift_amount = 0;
12044c1414c8SBarry Smith   sctx.nshift       = 0;
12054c1414c8SBarry Smith 
12064c1414c8SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
12074c1414c8SBarry Smith   ierr  = ISGetIndices(iscol,&c);CHKERRQ(ierr);
12084c1414c8SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
12094c1414c8SBarry Smith   ierr  = PetscMalloc((3*n+1)*sizeof(PetscScalar),&rtmp1);CHKERRQ(ierr);
12104c1414c8SBarry Smith   ierr  = PetscMemzero(rtmp1,(3*n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
12114c1414c8SBarry Smith   ics   = ic ;
12124c1414c8SBarry Smith   rtmp2 = rtmp1 + n;
12134c1414c8SBarry Smith   rtmp3 = rtmp2 + n;
12144c1414c8SBarry Smith 
12154c1414c8SBarry Smith   node_max = a->inode.node_count;
12164c1414c8SBarry Smith   ns       = a->inode.size ;
12174c1414c8SBarry Smith   if (!ns){
12184c1414c8SBarry Smith     SETERRQ(PETSC_ERR_PLIB,"Matrix without inode information");
12194c1414c8SBarry Smith   }
12204c1414c8SBarry Smith 
12214c1414c8SBarry Smith   /* If max inode size > 3, split it into two inodes.*/
12224c1414c8SBarry Smith   /* also map the inode sizes according to the ordering */
12234c1414c8SBarry Smith   ierr = PetscMalloc((n+1)* sizeof(PetscInt),&tmp_vec1);CHKERRQ(ierr);
12244c1414c8SBarry Smith   for (i=0,j=0; i<node_max; ++i,++j){
12254c1414c8SBarry Smith     if (ns[i]>3) {
12264c1414c8SBarry Smith       tmp_vec1[j] = ns[i]/2; /* Assuming ns[i] < =5  */
12274c1414c8SBarry Smith       ++j;
12284c1414c8SBarry Smith       tmp_vec1[j] = ns[i] - tmp_vec1[j-1];
12294c1414c8SBarry Smith     } else {
12304c1414c8SBarry Smith       tmp_vec1[j] = ns[i];
12314c1414c8SBarry Smith     }
12324c1414c8SBarry Smith   }
12334c1414c8SBarry Smith   /* Use the correct node_max */
12344c1414c8SBarry Smith   node_max = j;
12354c1414c8SBarry Smith 
12364c1414c8SBarry Smith   /* Now reorder the inode info based on mat re-ordering info */
12374c1414c8SBarry Smith   /* First create a row -> inode_size_array_index map */
12384c1414c8SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt)+1,&nsmap);CHKERRQ(ierr);
12394c1414c8SBarry Smith   ierr = PetscMalloc(node_max*sizeof(PetscInt)+1,&tmp_vec2);CHKERRQ(ierr);
12404c1414c8SBarry Smith   for (i=0,row=0; i<node_max; i++) {
12414c1414c8SBarry Smith     nodesz = tmp_vec1[i];
12424c1414c8SBarry Smith     for (j=0; j<nodesz; j++,row++) {
12434c1414c8SBarry Smith       nsmap[row] = i;
12444c1414c8SBarry Smith     }
12454c1414c8SBarry Smith   }
12464c1414c8SBarry Smith   /* Using nsmap, create a reordered ns structure */
12474c1414c8SBarry Smith   for (i=0,j=0; i< node_max; i++) {
12484c1414c8SBarry Smith     nodesz       = tmp_vec1[nsmap[r[j]]];    /* here the reordered row_no is in r[] */
12494c1414c8SBarry Smith     tmp_vec2[i]  = nodesz;
12504c1414c8SBarry Smith     j           += nodesz;
12514c1414c8SBarry Smith   }
12524c1414c8SBarry Smith   ierr = PetscFree(nsmap);CHKERRQ(ierr);
12534c1414c8SBarry Smith   ierr = PetscFree(tmp_vec1);CHKERRQ(ierr);
12544c1414c8SBarry Smith   /* Now use the correct ns */
12554c1414c8SBarry Smith   ns = tmp_vec2;
12564c1414c8SBarry Smith 
12574c1414c8SBarry Smith   do {
12584c1414c8SBarry Smith     sctx.lushift = PETSC_FALSE;
12594c1414c8SBarry Smith     /* Now loop over each block-row, and do the factorization */
12604c1414c8SBarry Smith     for (i=0,row=0; i<node_max; i++) {
12614c1414c8SBarry Smith       nodesz = ns[i];
12624c1414c8SBarry Smith       nz     = bi[row+1] - bi[row];
12634c1414c8SBarry Smith       bjtmp  = bj + bi[row];
12644c1414c8SBarry Smith 
12654c1414c8SBarry Smith       switch (nodesz){
12664c1414c8SBarry Smith       case 1:
12674c1414c8SBarry Smith         for  (j=0; j<nz; j++){
12684c1414c8SBarry Smith           idx        = bjtmp[j];
12694c1414c8SBarry Smith           rtmp1[idx] = 0.0;
12704c1414c8SBarry Smith         }
12714c1414c8SBarry Smith 
12724c1414c8SBarry Smith         /* load in initial (unfactored row) */
12734c1414c8SBarry Smith         idx    = r[row];
12744c1414c8SBarry Smith         nz_tmp = ai[idx+1] - ai[idx];
12754c1414c8SBarry Smith         ajtmp  = aj + ai[idx];
12764c1414c8SBarry Smith         v1     = aa + ai[idx];
12774c1414c8SBarry Smith 
12784c1414c8SBarry Smith         for (j=0; j<nz_tmp; j++) {
12794c1414c8SBarry Smith           idx        = ics[ajtmp[j]];
12804c1414c8SBarry Smith           rtmp1[idx] = v1[j];
12814c1414c8SBarry Smith         }
12824c1414c8SBarry Smith         rtmp1[ics[r[row]]] += sctx.shift_amount;
12834c1414c8SBarry Smith 
12844c1414c8SBarry Smith         prow = *bjtmp++ ;
12854c1414c8SBarry Smith         while (prow < row) {
12864c1414c8SBarry Smith           pc1 = rtmp1 + prow;
12874c1414c8SBarry Smith           if (*pc1 != 0.0){
12884c1414c8SBarry Smith             pv   = ba + bd[prow];
12894c1414c8SBarry Smith             pj   = nbj + bd[prow];
12904c1414c8SBarry Smith             mul1 = *pc1 * *pv++;
12914c1414c8SBarry Smith             *pc1 = mul1;
12924c1414c8SBarry Smith             nz_tmp = bi[prow+1] - bd[prow] - 1;
12934c1414c8SBarry Smith             ierr = PetscLogFlops(2*nz_tmp);CHKERRQ(ierr);
12944c1414c8SBarry Smith             for (j=0; j<nz_tmp; j++) {
12954c1414c8SBarry Smith               tmp = pv[j];
12964c1414c8SBarry Smith               idx = pj[j];
12974c1414c8SBarry Smith               rtmp1[idx] -= mul1 * tmp;
12984c1414c8SBarry Smith             }
12994c1414c8SBarry Smith           }
13004c1414c8SBarry Smith           prow = *bjtmp++ ;
13014c1414c8SBarry Smith         }
13024c1414c8SBarry Smith         pj  = bj + bi[row];
13034c1414c8SBarry Smith         pc1 = ba + bi[row];
13044c1414c8SBarry Smith 
13054c1414c8SBarry Smith         sctx.pv    = rtmp1[row];
13064c1414c8SBarry Smith         rtmp1[row] = 1.0/rtmp1[row]; /* invert diag */
13074c1414c8SBarry Smith         rs         = 0.0;
13084c1414c8SBarry Smith         for (j=0; j<nz; j++) {
13094c1414c8SBarry Smith           idx    = pj[j];
13104c1414c8SBarry Smith           pc1[j] = rtmp1[idx]; /* rtmp1 -> ba */
13114c1414c8SBarry Smith           if (idx != row) rs += PetscAbsScalar(pc1[j]);
13124c1414c8SBarry Smith         }
13134c1414c8SBarry Smith         sctx.rs  = rs;
1314*c6b1f410SHong Zhang         ierr = MatLUCheckShift_inline(info,sctx,row,newshift);CHKERRQ(ierr);
13154c1414c8SBarry Smith         if (newshift == 1) goto endofwhile;
13164c1414c8SBarry Smith         break;
13174c1414c8SBarry Smith 
13184c1414c8SBarry Smith       case 2:
13194c1414c8SBarry Smith         for (j=0; j<nz; j++) {
13204c1414c8SBarry Smith           idx        = bjtmp[j];
13214c1414c8SBarry Smith           rtmp1[idx] = 0.0;
13224c1414c8SBarry Smith           rtmp2[idx] = 0.0;
13234c1414c8SBarry Smith         }
13244c1414c8SBarry Smith 
13254c1414c8SBarry Smith         /* load in initial (unfactored row) */
13264c1414c8SBarry Smith         idx    = r[row];
13274c1414c8SBarry Smith         nz_tmp = ai[idx+1] - ai[idx];
13284c1414c8SBarry Smith         ajtmp  = aj + ai[idx];
13294c1414c8SBarry Smith         v1     = aa + ai[idx];
13304c1414c8SBarry Smith         v2     = aa + ai[idx+1];
13314c1414c8SBarry Smith         for (j=0; j<nz_tmp; j++) {
13324c1414c8SBarry Smith           idx        = ics[ajtmp[j]];
13334c1414c8SBarry Smith           rtmp1[idx] = v1[j];
13344c1414c8SBarry Smith           rtmp2[idx] = v2[j];
13354c1414c8SBarry Smith         }
13364c1414c8SBarry Smith         rtmp1[ics[r[row]]]   += sctx.shift_amount;
13374c1414c8SBarry Smith         rtmp2[ics[r[row+1]]] += sctx.shift_amount;
13384c1414c8SBarry Smith 
13394c1414c8SBarry Smith         prow = *bjtmp++ ;
13404c1414c8SBarry Smith         while (prow < row) {
13414c1414c8SBarry Smith           pc1 = rtmp1 + prow;
13424c1414c8SBarry Smith           pc2 = rtmp2 + prow;
13434c1414c8SBarry Smith           if (*pc1 != 0.0 || *pc2 != 0.0){
13444c1414c8SBarry Smith             pv   = ba + bd[prow];
13454c1414c8SBarry Smith             pj   = nbj + bd[prow];
13464c1414c8SBarry Smith             mul1 = *pc1 * *pv;
13474c1414c8SBarry Smith             mul2 = *pc2 * *pv;
13484c1414c8SBarry Smith             ++pv;
13494c1414c8SBarry Smith             *pc1 = mul1;
13504c1414c8SBarry Smith             *pc2 = mul2;
13514c1414c8SBarry Smith 
13524c1414c8SBarry Smith             nz_tmp = bi[prow+1] - bd[prow] - 1;
13534c1414c8SBarry Smith             for (j=0; j<nz_tmp; j++) {
13544c1414c8SBarry Smith               tmp = pv[j];
13554c1414c8SBarry Smith               idx = pj[j];
13564c1414c8SBarry Smith               rtmp1[idx] -= mul1 * tmp;
13574c1414c8SBarry Smith               rtmp2[idx] -= mul2 * tmp;
13584c1414c8SBarry Smith             }
13594c1414c8SBarry Smith             ierr = PetscLogFlops(4*nz_tmp);CHKERRQ(ierr);
13604c1414c8SBarry Smith           }
13614c1414c8SBarry Smith           prow = *bjtmp++ ;
13624c1414c8SBarry Smith         }
13634c1414c8SBarry Smith 
13644c1414c8SBarry Smith         /* Now take care of diagonal 2x2 block. Note: prow = row here */
13654c1414c8SBarry Smith         pc1 = rtmp1 + prow;
13664c1414c8SBarry Smith         pc2 = rtmp2 + prow;
13674c1414c8SBarry Smith 
13684c1414c8SBarry Smith         sctx.pv = *pc1;
13694c1414c8SBarry Smith         pj      = bj + bi[prow];
13704c1414c8SBarry Smith         rs      = 0.0;
13714c1414c8SBarry Smith         for (j=0; j<nz; j++){
13724c1414c8SBarry Smith           idx = pj[j];
13734c1414c8SBarry Smith           if (idx != prow) rs += PetscAbsScalar(rtmp1[idx]);
13744c1414c8SBarry Smith         }
13754c1414c8SBarry Smith         sctx.rs = rs;
1376*c6b1f410SHong Zhang         ierr = MatLUCheckShift_inline(info,sctx,row,newshift);CHKERRQ(ierr);
13774c1414c8SBarry Smith         if (newshift == 1) goto endofwhile;
13784c1414c8SBarry Smith 
13794c1414c8SBarry Smith         if (*pc2 != 0.0){
13804c1414c8SBarry Smith           pj     = nbj + bd[prow];
13814c1414c8SBarry Smith           mul2   = (*pc2)/(*pc1); /* since diag is not yet inverted.*/
13824c1414c8SBarry Smith           *pc2   = mul2;
13834c1414c8SBarry Smith           nz_tmp = bi[prow+1] - bd[prow] - 1;
13844c1414c8SBarry Smith           for (j=0; j<nz_tmp; j++) {
13854c1414c8SBarry Smith             idx = pj[j] ;
13864c1414c8SBarry Smith             tmp = rtmp1[idx];
13874c1414c8SBarry Smith             rtmp2[idx] -= mul2 * tmp;
13884c1414c8SBarry Smith           }
13894c1414c8SBarry Smith           ierr = PetscLogFlops(2*nz_tmp);CHKERRQ(ierr);
13904c1414c8SBarry Smith         }
13914c1414c8SBarry Smith 
13924c1414c8SBarry Smith         pj  = bj + bi[row];
13934c1414c8SBarry Smith         pc1 = ba + bi[row];
13944c1414c8SBarry Smith         pc2 = ba + bi[row+1];
13954c1414c8SBarry Smith 
13964c1414c8SBarry Smith         sctx.pv = rtmp2[row+1];
13974c1414c8SBarry Smith         rs = 0.0;
13984c1414c8SBarry Smith         rtmp1[row]   = 1.0/rtmp1[row];
13994c1414c8SBarry Smith         rtmp2[row+1] = 1.0/rtmp2[row+1];
14004c1414c8SBarry Smith         /* copy row entries from dense representation to sparse */
14014c1414c8SBarry Smith         for (j=0; j<nz; j++) {
14024c1414c8SBarry Smith           idx    = pj[j];
14034c1414c8SBarry Smith           pc1[j] = rtmp1[idx];
14044c1414c8SBarry Smith           pc2[j] = rtmp2[idx];
14054c1414c8SBarry Smith           if (idx != row+1) rs += PetscAbsScalar(pc2[j]);
14064c1414c8SBarry Smith         }
14074c1414c8SBarry Smith         sctx.rs = rs;
1408*c6b1f410SHong Zhang         ierr = MatLUCheckShift_inline(info,sctx,row+1,newshift);CHKERRQ(ierr);
14094c1414c8SBarry Smith         if (newshift == 1) goto endofwhile;
14104c1414c8SBarry Smith         break;
14114c1414c8SBarry Smith 
14124c1414c8SBarry Smith       case 3:
14134c1414c8SBarry Smith         for  (j=0; j<nz; j++) {
14144c1414c8SBarry Smith           idx        = bjtmp[j];
14154c1414c8SBarry Smith           rtmp1[idx] = 0.0;
14164c1414c8SBarry Smith           rtmp2[idx] = 0.0;
14174c1414c8SBarry Smith           rtmp3[idx] = 0.0;
14184c1414c8SBarry Smith         }
14194c1414c8SBarry Smith         /* copy the nonzeros for the 3 rows from sparse representation to dense in rtmp*[] */
14204c1414c8SBarry Smith         idx    = r[row];
14214c1414c8SBarry Smith         nz_tmp = ai[idx+1] - ai[idx];
14224c1414c8SBarry Smith         ajtmp = aj + ai[idx];
14234c1414c8SBarry Smith         v1    = aa + ai[idx];
14244c1414c8SBarry Smith         v2    = aa + ai[idx+1];
14254c1414c8SBarry Smith         v3    = aa + ai[idx+2];
14264c1414c8SBarry Smith         for (j=0; j<nz_tmp; j++) {
14274c1414c8SBarry Smith           idx        = ics[ajtmp[j]];
14284c1414c8SBarry Smith           rtmp1[idx] = v1[j];
14294c1414c8SBarry Smith           rtmp2[idx] = v2[j];
14304c1414c8SBarry Smith           rtmp3[idx] = v3[j];
14314c1414c8SBarry Smith         }
14324c1414c8SBarry Smith         rtmp1[ics[r[row]]]   += sctx.shift_amount;
14334c1414c8SBarry Smith         rtmp2[ics[r[row+1]]] += sctx.shift_amount;
14344c1414c8SBarry Smith         rtmp3[ics[r[row+2]]] += sctx.shift_amount;
14354c1414c8SBarry Smith 
14364c1414c8SBarry Smith         /* loop over all pivot row blocks above this row block */
14374c1414c8SBarry Smith         prow = *bjtmp++ ;
14384c1414c8SBarry Smith         while (prow < row) {
14394c1414c8SBarry Smith           pc1 = rtmp1 + prow;
14404c1414c8SBarry Smith           pc2 = rtmp2 + prow;
14414c1414c8SBarry Smith           pc3 = rtmp3 + prow;
14424c1414c8SBarry Smith           if (*pc1 != 0.0 || *pc2 != 0.0 || *pc3 !=0.0){
14434c1414c8SBarry Smith             pv   = ba  + bd[prow];
14444c1414c8SBarry Smith             pj   = nbj + bd[prow];
14454c1414c8SBarry Smith             mul1 = *pc1 * *pv;
14464c1414c8SBarry Smith             mul2 = *pc2 * *pv;
14474c1414c8SBarry Smith             mul3 = *pc3 * *pv;
14484c1414c8SBarry Smith             ++pv;
14494c1414c8SBarry Smith             *pc1 = mul1;
14504c1414c8SBarry Smith             *pc2 = mul2;
14514c1414c8SBarry Smith             *pc3 = mul3;
14524c1414c8SBarry Smith 
14534c1414c8SBarry Smith             nz_tmp = bi[prow+1] - bd[prow] - 1;
14544c1414c8SBarry Smith             /* update this row based on pivot row */
14554c1414c8SBarry Smith             for (j=0; j<nz_tmp; j++) {
14564c1414c8SBarry Smith               tmp = pv[j];
14574c1414c8SBarry Smith               idx = pj[j];
14584c1414c8SBarry Smith               rtmp1[idx] -= mul1 * tmp;
14594c1414c8SBarry Smith               rtmp2[idx] -= mul2 * tmp;
14604c1414c8SBarry Smith               rtmp3[idx] -= mul3 * tmp;
14614c1414c8SBarry Smith             }
14624c1414c8SBarry Smith             ierr = PetscLogFlops(6*nz_tmp);CHKERRQ(ierr);
14634c1414c8SBarry Smith           }
14644c1414c8SBarry Smith           prow = *bjtmp++ ;
14654c1414c8SBarry Smith         }
14664c1414c8SBarry Smith 
14674c1414c8SBarry Smith         /* Now take care of diagonal 3x3 block in this set of rows */
14684c1414c8SBarry Smith         /* note: prow = row here */
14694c1414c8SBarry Smith         pc1 = rtmp1 + prow;
14704c1414c8SBarry Smith         pc2 = rtmp2 + prow;
14714c1414c8SBarry Smith         pc3 = rtmp3 + prow;
14724c1414c8SBarry Smith 
14734c1414c8SBarry Smith         sctx.pv = *pc1;
14744c1414c8SBarry Smith         pj      = bj + bi[prow];
14754c1414c8SBarry Smith         rs      = 0.0;
14764c1414c8SBarry Smith         for (j=0; j<nz; j++){
14774c1414c8SBarry Smith           idx = pj[j];
14784c1414c8SBarry Smith           if (idx != row) rs += PetscAbsScalar(rtmp1[idx]);
14794c1414c8SBarry Smith         }
14804c1414c8SBarry Smith         sctx.rs = rs;
1481*c6b1f410SHong Zhang         ierr = MatLUCheckShift_inline(info,sctx,row,newshift);CHKERRQ(ierr);
14824c1414c8SBarry Smith         if (newshift == 1) goto endofwhile;
14834c1414c8SBarry Smith 
14844c1414c8SBarry Smith         if (*pc2 != 0.0 || *pc3 != 0.0){
14854c1414c8SBarry Smith           mul2 = (*pc2)/(*pc1);
14864c1414c8SBarry Smith           mul3 = (*pc3)/(*pc1);
14874c1414c8SBarry Smith           *pc2 = mul2;
14884c1414c8SBarry Smith           *pc3 = mul3;
14894c1414c8SBarry Smith           nz_tmp = bi[prow+1] - bd[prow] - 1;
14904c1414c8SBarry Smith           pj     = nbj + bd[prow];
14914c1414c8SBarry Smith           for (j=0; j<nz_tmp; j++) {
14924c1414c8SBarry Smith             idx = pj[j] ;
14934c1414c8SBarry Smith             tmp = rtmp1[idx];
14944c1414c8SBarry Smith             rtmp2[idx] -= mul2 * tmp;
14954c1414c8SBarry Smith             rtmp3[idx] -= mul3 * tmp;
14964c1414c8SBarry Smith           }
14974c1414c8SBarry Smith           ierr = PetscLogFlops(4*nz_tmp);CHKERRQ(ierr);
14984c1414c8SBarry Smith         }
14994c1414c8SBarry Smith         ++prow;
15004c1414c8SBarry Smith 
15014c1414c8SBarry Smith         pc2 = rtmp2 + prow;
15024c1414c8SBarry Smith         pc3 = rtmp3 + prow;
15034c1414c8SBarry Smith         sctx.pv = *pc2;
15044c1414c8SBarry Smith         pj      = bj + bi[prow];
15054c1414c8SBarry Smith         rs      = 0.0;
15064c1414c8SBarry Smith         for (j=0; j<nz; j++){
15074c1414c8SBarry Smith           idx = pj[j];
15084c1414c8SBarry Smith           if (idx != prow) rs += PetscAbsScalar(rtmp2[idx]);
15094c1414c8SBarry Smith         }
15104c1414c8SBarry Smith         sctx.rs = rs;
1511*c6b1f410SHong Zhang         ierr = MatLUCheckShift_inline(info,sctx,row+1,newshift);CHKERRQ(ierr);
15124c1414c8SBarry Smith         if (newshift == 1) goto endofwhile;
15134c1414c8SBarry Smith 
15144c1414c8SBarry Smith         if (*pc3 != 0.0){
15154c1414c8SBarry Smith           mul3   = (*pc3)/(*pc2);
15164c1414c8SBarry Smith           *pc3   = mul3;
15174c1414c8SBarry Smith           pj     = nbj + bd[prow];
15184c1414c8SBarry Smith           nz_tmp = bi[prow+1] - bd[prow] - 1;
15194c1414c8SBarry Smith           for (j=0; j<nz_tmp; j++) {
15204c1414c8SBarry Smith             idx = pj[j] ;
15214c1414c8SBarry Smith             tmp = rtmp2[idx];
15224c1414c8SBarry Smith             rtmp3[idx] -= mul3 * tmp;
15234c1414c8SBarry Smith           }
15244c1414c8SBarry Smith           ierr = PetscLogFlops(4*nz_tmp);CHKERRQ(ierr);
15254c1414c8SBarry Smith         }
15264c1414c8SBarry Smith 
15274c1414c8SBarry Smith         pj  = bj + bi[row];
15284c1414c8SBarry Smith         pc1 = ba + bi[row];
15294c1414c8SBarry Smith         pc2 = ba + bi[row+1];
15304c1414c8SBarry Smith         pc3 = ba + bi[row+2];
15314c1414c8SBarry Smith 
15324c1414c8SBarry Smith         sctx.pv = rtmp3[row+2];
15334c1414c8SBarry Smith         rs = 0.0;
15344c1414c8SBarry Smith         rtmp1[row]   = 1.0/rtmp1[row];
15354c1414c8SBarry Smith         rtmp2[row+1] = 1.0/rtmp2[row+1];
15364c1414c8SBarry Smith         rtmp3[row+2] = 1.0/rtmp3[row+2];
15374c1414c8SBarry Smith         /* copy row entries from dense representation to sparse */
15384c1414c8SBarry Smith         for (j=0; j<nz; j++) {
15394c1414c8SBarry Smith           idx    = pj[j];
15404c1414c8SBarry Smith           pc1[j] = rtmp1[idx];
15414c1414c8SBarry Smith           pc2[j] = rtmp2[idx];
15424c1414c8SBarry Smith           pc3[j] = rtmp3[idx];
15434c1414c8SBarry Smith           if (idx != row+2) rs += PetscAbsScalar(pc3[j]);
15444c1414c8SBarry Smith         }
15454c1414c8SBarry Smith 
15464c1414c8SBarry Smith         sctx.rs = rs;
1547*c6b1f410SHong Zhang         ierr = MatLUCheckShift_inline(info,sctx,row+2,newshift);CHKERRQ(ierr);
15484c1414c8SBarry Smith         if (newshift == 1) goto endofwhile;
15494c1414c8SBarry Smith         break;
15504c1414c8SBarry Smith 
15514c1414c8SBarry Smith       default:
15524c1414c8SBarry Smith         SETERRQ(PETSC_ERR_COR,"Node size not yet supported \n");
15534c1414c8SBarry Smith       }
15544c1414c8SBarry Smith       row += nodesz;                 /* Update the row */
15554c1414c8SBarry Smith     }
15564c1414c8SBarry Smith     endofwhile:;
15574c1414c8SBarry Smith   } while (sctx.lushift);
15584c1414c8SBarry Smith   ierr = PetscFree(rtmp1);CHKERRQ(ierr);
15594c1414c8SBarry Smith   ierr = PetscFree(tmp_vec2);CHKERRQ(ierr);
15604c1414c8SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
15614c1414c8SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
15624c1414c8SBarry Smith   ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
15634c1414c8SBarry Smith   C->factor      = FACTOR_LU;
15644c1414c8SBarry Smith   C->assembled   = PETSC_TRUE;
15654c1414c8SBarry Smith   if (sctx.nshift) {
15664c1414c8SBarry Smith     if (info->shiftnz) {
15674c1414c8SBarry Smith       ierr = PetscInfo2(0,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
15684c1414c8SBarry Smith     } else if (info->shiftpd) {
15694c1414c8SBarry Smith       ierr = PetscInfo4(0,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
15704c1414c8SBarry Smith     }
15714c1414c8SBarry Smith   }
15724c1414c8SBarry Smith   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
15734c1414c8SBarry Smith   PetscFunctionReturn(0);
15744c1414c8SBarry Smith }
15754c1414c8SBarry Smith 
15764c1414c8SBarry Smith /*
15774c1414c8SBarry Smith      Makes a longer coloring[] array and calls the usual code with that
15784c1414c8SBarry Smith */
15794c1414c8SBarry Smith #undef __FUNCT__
15804c1414c8SBarry Smith #define __FUNCT__ "MatColoringPatch_Inode"
15814c1414c8SBarry Smith PetscErrorCode MatColoringPatch_Inode(Mat mat,PetscInt ncolors,PetscInt nin,ISColoringValue coloring[],ISColoring *iscoloring)
15824c1414c8SBarry Smith {
15834c1414c8SBarry Smith   Mat_SeqAIJ       *a = (Mat_SeqAIJ*)mat->data;
15844c1414c8SBarry Smith   PetscErrorCode  ierr;
15854c1414c8SBarry Smith   PetscInt        n = mat->cmap.n,m = a->inode.node_count,j,*ns = a->inode.size,row;
15864c1414c8SBarry Smith   PetscInt        *colorused,i;
15874c1414c8SBarry Smith   ISColoringValue *newcolor;
15884c1414c8SBarry Smith 
15894c1414c8SBarry Smith   PetscFunctionBegin;
15904c1414c8SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&newcolor);CHKERRQ(ierr);
15914c1414c8SBarry Smith   /* loop over inodes, marking a color for each column*/
15924c1414c8SBarry Smith   row = 0;
15934c1414c8SBarry Smith   for (i=0; i<m; i++){
15944c1414c8SBarry Smith     for (j=0; j<ns[i]; j++) {
15954c1414c8SBarry Smith       newcolor[row++] = coloring[i] + j*ncolors;
15964c1414c8SBarry Smith     }
15974c1414c8SBarry Smith   }
15984c1414c8SBarry Smith 
15994c1414c8SBarry Smith   /* eliminate unneeded colors */
16004c1414c8SBarry Smith   ierr = PetscMalloc(5*ncolors*sizeof(PetscInt),&colorused);CHKERRQ(ierr);
16014c1414c8SBarry Smith   ierr = PetscMemzero(colorused,5*ncolors*sizeof(PetscInt));CHKERRQ(ierr);
16024c1414c8SBarry Smith   for (i=0; i<n; i++) {
16034c1414c8SBarry Smith     colorused[newcolor[i]] = 1;
16044c1414c8SBarry Smith   }
16054c1414c8SBarry Smith 
16064c1414c8SBarry Smith   for (i=1; i<5*ncolors; i++) {
16074c1414c8SBarry Smith     colorused[i] += colorused[i-1];
16084c1414c8SBarry Smith   }
16094c1414c8SBarry Smith   ncolors = colorused[5*ncolors-1];
16104c1414c8SBarry Smith   for (i=0; i<n; i++) {
16114c1414c8SBarry Smith     newcolor[i] = colorused[newcolor[i]];
16124c1414c8SBarry Smith   }
16134c1414c8SBarry Smith   ierr = PetscFree(colorused);CHKERRQ(ierr);
16144c1414c8SBarry Smith   ierr = ISColoringCreate(mat->comm,ncolors,n,newcolor,iscoloring);CHKERRQ(ierr);
16154c1414c8SBarry Smith   ierr = PetscFree(coloring);CHKERRQ(ierr);
16164c1414c8SBarry Smith   PetscFunctionReturn(0);
16174c1414c8SBarry Smith }
16184c1414c8SBarry Smith 
16194c1414c8SBarry Smith /*
16204c1414c8SBarry Smith     samestructure indicates that the matrix has not changed its nonzero structure so we
16214c1414c8SBarry Smith     do not need to recompute the inodes
16224c1414c8SBarry Smith */
16234c1414c8SBarry Smith #undef __FUNCT__
16244c1414c8SBarry Smith #define __FUNCT__ "Mat_CheckInode"
16254c1414c8SBarry Smith PetscErrorCode Mat_CheckInode(Mat A,PetscTruth samestructure)
16264c1414c8SBarry Smith {
16274c1414c8SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
16284c1414c8SBarry Smith   PetscErrorCode ierr;
16294c1414c8SBarry Smith   PetscInt       i,j,m,nzx,nzy,*idx,*idy,*ns,*ii,node_count,blk_size;
1630d74fe821SBarry Smith   PetscTruth     flag;
16314c1414c8SBarry Smith 
16324c1414c8SBarry Smith   PetscFunctionBegin;
1633d74fe821SBarry Smith   if (!a->inode.use)                     PetscFunctionReturn(0);
16344c1414c8SBarry Smith   if (a->inode.checked && samestructure) PetscFunctionReturn(0);
16354c1414c8SBarry Smith 
16364c1414c8SBarry Smith 
16374c1414c8SBarry Smith   m = A->rmap.n;
16384c1414c8SBarry Smith   if (a->inode.size) {ns = a->inode.size;}
16394c1414c8SBarry Smith   else {ierr = PetscMalloc((m+1)*sizeof(PetscInt),&ns);CHKERRQ(ierr);}
16404c1414c8SBarry Smith 
16414c1414c8SBarry Smith   i          = 0;
16424c1414c8SBarry Smith   node_count = 0;
16434c1414c8SBarry Smith   idx        = a->j;
16444c1414c8SBarry Smith   ii         = a->i;
16454c1414c8SBarry Smith   while (i < m){                /* For each row */
16464c1414c8SBarry Smith     nzx = ii[i+1] - ii[i];       /* Number of nonzeros */
16474c1414c8SBarry Smith     /* Limits the number of elements in a node to 'a->inode.limit' */
16484c1414c8SBarry Smith     for (j=i+1,idy=idx,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
16494c1414c8SBarry Smith       nzy     = ii[j+1] - ii[j]; /* Same number of nonzeros */
16504c1414c8SBarry Smith       if (nzy != nzx) break;
16514c1414c8SBarry Smith       idy  += nzx;             /* Same nonzero pattern */
16524c1414c8SBarry Smith       ierr = PetscMemcmp(idx,idy,nzx*sizeof(PetscInt),&flag);CHKERRQ(ierr);
16534c1414c8SBarry Smith       if (!flag) break;
16544c1414c8SBarry Smith     }
16554c1414c8SBarry Smith     ns[node_count++] = blk_size;
16564c1414c8SBarry Smith     idx += blk_size*nzx;
16574c1414c8SBarry Smith     i    = j;
16584c1414c8SBarry Smith   }
16594c1414c8SBarry Smith   /* If not enough inodes found,, do not use inode version of the routines */
16604c1414c8SBarry Smith   if (!a->inode.size && m && node_count > 0.9*m) {
16614c1414c8SBarry Smith     ierr = PetscFree(ns);CHKERRQ(ierr);
16624c1414c8SBarry Smith     a->inode.node_count     = 0;
16634c1414c8SBarry Smith     a->inode.size           = PETSC_NULL;
16644c1414c8SBarry Smith     a->inode.use            = PETSC_FALSE;
16654c1414c8SBarry Smith     ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr);
16664c1414c8SBarry Smith   } else {
16674c1414c8SBarry Smith     A->ops->mult            = MatMult_Inode;
16684c1414c8SBarry Smith     A->ops->multadd         = MatMultAdd_Inode;
16694c1414c8SBarry Smith     A->ops->solve           = MatSolve_Inode;
16704c1414c8SBarry Smith     A->ops->lufactornumeric = MatLUFactorNumeric_Inode;
16714c1414c8SBarry Smith     A->ops->getrowij        = MatGetRowIJ_Inode;
16724c1414c8SBarry Smith     A->ops->restorerowij    = MatRestoreRowIJ_Inode;
16734c1414c8SBarry Smith     A->ops->getcolumnij     = MatGetColumnIJ_Inode;
16744c1414c8SBarry Smith     A->ops->restorecolumnij = MatRestoreColumnIJ_Inode;
16754c1414c8SBarry Smith     A->ops->coloringpatch   = MatColoringPatch_Inode;
16764c1414c8SBarry Smith     a->inode.node_count     = node_count;
16774c1414c8SBarry Smith     a->inode.size           = ns;
16784c1414c8SBarry Smith     ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr);
16794c1414c8SBarry Smith   }
16804c1414c8SBarry Smith   PetscFunctionReturn(0);
16814c1414c8SBarry Smith }
16824c1414c8SBarry Smith 
16834c1414c8SBarry Smith /*
16844c1414c8SBarry Smith      This is really ugly. if inodes are used this replaces the
16854c1414c8SBarry Smith   permutations with ones that correspond to rows/cols of the matrix
16864c1414c8SBarry Smith   rather then inode blocks
16874c1414c8SBarry Smith */
16884c1414c8SBarry Smith #undef __FUNCT__
16894c1414c8SBarry Smith #define __FUNCT__ "MatInodeAdjustForInodes"
16904c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat A,IS *rperm,IS *cperm)
16914c1414c8SBarry Smith {
16924c1414c8SBarry Smith   PetscErrorCode ierr,(*f)(Mat,IS*,IS*);
16934c1414c8SBarry Smith 
16944c1414c8SBarry Smith   PetscFunctionBegin;
16954c1414c8SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatInodeAdjustForInodes_C",(void (**)(void))&f);CHKERRQ(ierr);
16964c1414c8SBarry Smith   if (f) {
16974c1414c8SBarry Smith     ierr = (*f)(A,rperm,cperm);CHKERRQ(ierr);
16984c1414c8SBarry Smith   }
16994c1414c8SBarry Smith   PetscFunctionReturn(0);
17004c1414c8SBarry Smith }
17014c1414c8SBarry Smith 
17024c1414c8SBarry Smith EXTERN_C_BEGIN
17034c1414c8SBarry Smith #undef __FUNCT__
17044c1414c8SBarry Smith #define __FUNCT__ "MatAdjustForInodes_Inode"
17054c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes_Inode(Mat A,IS *rperm,IS *cperm)
17064c1414c8SBarry Smith {
17074c1414c8SBarry Smith   Mat_SeqAIJ      *a=(Mat_SeqAIJ*)A->data;
17084c1414c8SBarry Smith   PetscErrorCode ierr;
17094c1414c8SBarry Smith   PetscInt       m = A->rmap.n,n = A->cmap.n,i,j,*ridx,*cidx,nslim_row = a->inode.node_count;
17104c1414c8SBarry Smith   PetscInt       row,col,*permr,*permc,*ns_row =  a->inode.size,*tns,start_val,end_val,indx;
17114c1414c8SBarry Smith   PetscInt       nslim_col,*ns_col;
17124c1414c8SBarry Smith   IS             ris = *rperm,cis = *cperm;
17134c1414c8SBarry Smith 
17144c1414c8SBarry Smith   PetscFunctionBegin;
17154c1414c8SBarry Smith   if (!a->inode.size) PetscFunctionReturn(0); /* no inodes so return */
17164c1414c8SBarry Smith   if (a->inode.node_count == m) PetscFunctionReturn(0); /* all inodes are of size 1 */
17174c1414c8SBarry Smith 
17184c1414c8SBarry Smith   ierr  = Mat_CreateColInode(A,&nslim_col,&ns_col);CHKERRQ(ierr);
17194c1414c8SBarry Smith   ierr  = PetscMalloc((((nslim_row>nslim_col)?nslim_row:nslim_col)+1)*sizeof(PetscInt),&tns);CHKERRQ(ierr);
17204c1414c8SBarry Smith   ierr  = PetscMalloc((m+n+1)*sizeof(PetscInt),&permr);CHKERRQ(ierr);
17214c1414c8SBarry Smith   permc = permr + m;
17224c1414c8SBarry Smith 
17234c1414c8SBarry Smith   ierr  = ISGetIndices(ris,&ridx);CHKERRQ(ierr);
17244c1414c8SBarry Smith   ierr  = ISGetIndices(cis,&cidx);CHKERRQ(ierr);
17254c1414c8SBarry Smith 
17264c1414c8SBarry Smith   /* Form the inode structure for the rows of permuted matric using inv perm*/
17274c1414c8SBarry Smith   for (i=0,tns[0]=0; i<nslim_row; ++i) tns[i+1] = tns[i] + ns_row[i];
17284c1414c8SBarry Smith 
17294c1414c8SBarry Smith   /* Construct the permutations for rows*/
17304c1414c8SBarry Smith   for (i=0,row = 0; i<nslim_row; ++i){
17314c1414c8SBarry Smith     indx      = ridx[i];
17324c1414c8SBarry Smith     start_val = tns[indx];
17334c1414c8SBarry Smith     end_val   = tns[indx + 1];
17344c1414c8SBarry Smith     for (j=start_val; j<end_val; ++j,++row) permr[row]= j;
17354c1414c8SBarry Smith   }
17364c1414c8SBarry Smith 
17374c1414c8SBarry Smith   /* Form the inode structure for the columns of permuted matrix using inv perm*/
17384c1414c8SBarry Smith   for (i=0,tns[0]=0; i<nslim_col; ++i) tns[i+1] = tns[i] + ns_col[i];
17394c1414c8SBarry Smith 
17404c1414c8SBarry Smith  /* Construct permutations for columns */
17414c1414c8SBarry Smith   for (i=0,col=0; i<nslim_col; ++i){
17424c1414c8SBarry Smith     indx      = cidx[i];
17434c1414c8SBarry Smith     start_val = tns[indx];
17444c1414c8SBarry Smith     end_val   = tns[indx + 1];
17454c1414c8SBarry Smith     for (j = start_val; j<end_val; ++j,++col) permc[col]= j;
17464c1414c8SBarry Smith   }
17474c1414c8SBarry Smith 
17484c1414c8SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,permr,rperm);CHKERRQ(ierr);
17494c1414c8SBarry Smith   ISSetPermutation(*rperm);
17504c1414c8SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,permc,cperm);CHKERRQ(ierr);
17514c1414c8SBarry Smith   ISSetPermutation(*cperm);
17524c1414c8SBarry Smith 
17534c1414c8SBarry Smith   ierr  = ISRestoreIndices(ris,&ridx);CHKERRQ(ierr);
17544c1414c8SBarry Smith   ierr  = ISRestoreIndices(cis,&cidx);CHKERRQ(ierr);
17554c1414c8SBarry Smith 
17564c1414c8SBarry Smith   ierr = PetscFree(ns_col);CHKERRQ(ierr);
17574c1414c8SBarry Smith   ierr = PetscFree(permr);CHKERRQ(ierr);
17584c1414c8SBarry Smith   ierr = ISDestroy(cis);CHKERRQ(ierr);
17594c1414c8SBarry Smith   ierr = ISDestroy(ris);CHKERRQ(ierr);
17604c1414c8SBarry Smith   ierr = PetscFree(tns);CHKERRQ(ierr);
17614c1414c8SBarry Smith   PetscFunctionReturn(0);
17624c1414c8SBarry Smith }
17634c1414c8SBarry Smith EXTERN_C_END
17644c1414c8SBarry Smith 
17654c1414c8SBarry Smith #undef __FUNCT__
17664c1414c8SBarry Smith #define __FUNCT__ "MatInodeGetInodeSizes"
17674c1414c8SBarry Smith /*@C
17684c1414c8SBarry Smith    MatInodeGetInodeSizes - Returns the inode information of the Inode matrix.
17694c1414c8SBarry Smith 
17704c1414c8SBarry Smith    Collective on Mat
17714c1414c8SBarry Smith 
17724c1414c8SBarry Smith    Input Parameter:
17734c1414c8SBarry Smith .  A - the Inode matrix or matrix derived from the Inode class -- e.g., SeqAIJ
17744c1414c8SBarry Smith 
17754c1414c8SBarry Smith    Output Parameter:
17764c1414c8SBarry Smith +  node_count - no of inodes present in the matrix.
17774c1414c8SBarry Smith .  sizes      - an array of size node_count,with sizes of each inode.
17784c1414c8SBarry Smith -  limit      - the max size used to generate the inodes.
17794c1414c8SBarry Smith 
17804c1414c8SBarry Smith    Level: advanced
17814c1414c8SBarry Smith 
17824c1414c8SBarry Smith    Notes: This routine returns some internal storage information
17834c1414c8SBarry Smith    of the matrix, it is intended to be used by advanced users.
17844c1414c8SBarry Smith    It should be called after the matrix is assembled.
17854c1414c8SBarry Smith    The contents of the sizes[] array should not be changed.
17864c1414c8SBarry Smith    PETSC_NULL may be passed for information not requested.
17874c1414c8SBarry Smith 
17884c1414c8SBarry Smith .keywords: matrix, seqaij, get, inode
17894c1414c8SBarry Smith 
17904c1414c8SBarry Smith .seealso: MatGetInfo()
17914c1414c8SBarry Smith @*/
17924c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat A,PetscInt *node_count,PetscInt *sizes[],PetscInt *limit)
17934c1414c8SBarry Smith {
17944c1414c8SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*[],PetscInt*);
17954c1414c8SBarry Smith 
17964c1414c8SBarry Smith   PetscFunctionBegin;
17974c1414c8SBarry Smith   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
17984c1414c8SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatInodeGetInodeSizes_C",(void (**)(void))&f);CHKERRQ(ierr);
17994c1414c8SBarry Smith   if (f) {
18004c1414c8SBarry Smith     ierr = (*f)(A,node_count,sizes,limit);CHKERRQ(ierr);
18014c1414c8SBarry Smith   }
18024c1414c8SBarry Smith   PetscFunctionReturn(0);
18034c1414c8SBarry Smith }
18044c1414c8SBarry Smith 
18054c1414c8SBarry Smith EXTERN_C_BEGIN
18064c1414c8SBarry Smith #undef __FUNCT__
18074c1414c8SBarry Smith #define __FUNCT__ "MatInodeGetInodeSizes_Inode"
18084c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes_Inode(Mat A,PetscInt *node_count,PetscInt *sizes[],PetscInt *limit)
18094c1414c8SBarry Smith {
18104c1414c8SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
18114c1414c8SBarry Smith 
18124c1414c8SBarry Smith   PetscFunctionBegin;
18134c1414c8SBarry Smith   if (node_count) *node_count = a->inode.node_count;
18144c1414c8SBarry Smith   if (sizes)      *sizes      = a->inode.size;
18154c1414c8SBarry Smith   if (limit)      *limit      = a->inode.limit;
18164c1414c8SBarry Smith   PetscFunctionReturn(0);
18174c1414c8SBarry Smith }
18184c1414c8SBarry Smith EXTERN_C_END
1819