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; 18d0f46423SBarry Smith n = A->cmap->n; 19d0f46423SBarry 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; 71d0f46423SBarry Smith m = A->rmap->n; 72d0f46423SBarry 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; 163d0f46423SBarry 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" 2348f7157efSSatish Balay static PetscErrorCode MatGetRowIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,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); 2428f7157efSSatish Balay if (!blockcompressed) { 2438f7157efSSatish Balay ierr = MatGetRowIJ_SeqAIJ(A,oshift,symmetric,blockcompressed,n,ia,ja,done);CHKERRQ(ierr);; 2448f7157efSSatish Balay } else if (symmetric) { 2454c1414c8SBarry Smith ierr = MatGetRowIJ_Inode_Symmetric(A,ia,ja,0,oshift);CHKERRQ(ierr); 2464c1414c8SBarry Smith } else { 2474c1414c8SBarry Smith ierr = MatGetRowIJ_Inode_Nonsymmetric(A,ia,ja,0,oshift);CHKERRQ(ierr); 2484c1414c8SBarry Smith } 2494c1414c8SBarry Smith PetscFunctionReturn(0); 2504c1414c8SBarry Smith } 2514c1414c8SBarry Smith 2524c1414c8SBarry Smith #undef __FUNCT__ 2534c1414c8SBarry Smith #define __FUNCT__ "MatRestoreRowIJ_Inode" 2548f7157efSSatish Balay static PetscErrorCode MatRestoreRowIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 2554c1414c8SBarry Smith { 2564c1414c8SBarry Smith PetscErrorCode ierr; 2574c1414c8SBarry Smith 2584c1414c8SBarry Smith PetscFunctionBegin; 2594c1414c8SBarry Smith if (!ia) PetscFunctionReturn(0); 2608f7157efSSatish Balay 2618f7157efSSatish Balay if (!blockcompressed) { 2628f7157efSSatish Balay ierr = MatRestoreRowIJ_SeqAIJ(A,oshift,symmetric,blockcompressed,n,ia,ja,done);CHKERRQ(ierr);; 2638f7157efSSatish Balay } else { 2644c1414c8SBarry Smith ierr = PetscFree(*ia);CHKERRQ(ierr); 2654c1414c8SBarry Smith ierr = PetscFree(*ja);CHKERRQ(ierr); 2668f7157efSSatish Balay } 2678f7157efSSatish Balay 2684c1414c8SBarry Smith PetscFunctionReturn(0); 2694c1414c8SBarry Smith } 2704c1414c8SBarry Smith 2714c1414c8SBarry Smith /* ----------------------------------------------------------- */ 2724c1414c8SBarry Smith 2734c1414c8SBarry Smith #undef __FUNCT__ 2744c1414c8SBarry Smith #define __FUNCT__ "MatGetColumnIJ_Inode_Nonsymmetric" 2754c1414c8SBarry Smith static PetscErrorCode MatGetColumnIJ_Inode_Nonsymmetric(Mat A,PetscInt *iia[],PetscInt *jja[],PetscInt ishift,PetscInt oshift) 2764c1414c8SBarry Smith { 2774c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2784c1414c8SBarry Smith PetscErrorCode ierr; 2794c1414c8SBarry Smith PetscInt *work,*ia,*ja,*j,nz,nslim_row, n,row,col,*ns_col,nslim_col; 2804c1414c8SBarry Smith PetscInt *tns,*tvc,*ns_row = a->inode.size,nsz,i1,i2,*ai= a->i,*aj = a->j; 2814c1414c8SBarry Smith 2824c1414c8SBarry Smith PetscFunctionBegin; 2834c1414c8SBarry Smith nslim_row = a->inode.node_count; 284d0f46423SBarry Smith n = A->cmap->n; 2854c1414c8SBarry Smith 2864c1414c8SBarry Smith /* Create The column_inode for this matrix */ 2874c1414c8SBarry Smith ierr = Mat_CreateColInode(A,&nslim_col,&ns_col);CHKERRQ(ierr); 2884c1414c8SBarry Smith 2894c1414c8SBarry Smith /* allocate space for reformated column_inode structure */ 2904c1414c8SBarry Smith ierr = PetscMalloc((nslim_col + 1)*sizeof(PetscInt),&tns);CHKERRQ(ierr); 2914c1414c8SBarry Smith ierr = PetscMalloc((n + 1)*sizeof(PetscInt),&tvc);CHKERRQ(ierr); 2924c1414c8SBarry Smith for (i1=0,tns[0]=0; i1<nslim_col; ++i1) tns[i1+1] = tns[i1] + ns_col[i1]; 2934c1414c8SBarry Smith 2944c1414c8SBarry Smith for (i1=0,col=0; i1<nslim_col; ++i1){ 2954c1414c8SBarry Smith nsz = ns_col[i1]; 2964c1414c8SBarry Smith for (i2=0; i2<nsz; ++i2,++col) 2974c1414c8SBarry Smith tvc[col] = i1; 2984c1414c8SBarry Smith } 2994c1414c8SBarry Smith /* allocate space for column pointers */ 3004c1414c8SBarry Smith ierr = PetscMalloc((nslim_col+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 3014c1414c8SBarry Smith *iia = ia; 3024c1414c8SBarry Smith ierr = PetscMemzero(ia,(nslim_col+1)*sizeof(PetscInt));CHKERRQ(ierr); 3034c1414c8SBarry Smith ierr = PetscMalloc((nslim_col+1)*sizeof(PetscInt),&work);CHKERRQ(ierr); 3044c1414c8SBarry Smith 3054c1414c8SBarry Smith /* determine the number of columns in each row */ 3064c1414c8SBarry Smith ia[0] = oshift; 3074c1414c8SBarry Smith for (i1=0,row=0; i1<nslim_row; row+=ns_row[i1],i1++) { 3084c1414c8SBarry Smith j = aj + ai[row] + ishift; 3094c1414c8SBarry Smith col = *j++ + ishift; 3104c1414c8SBarry Smith i2 = tvc[col]; 3114c1414c8SBarry Smith nz = ai[row+1] - ai[row]; 3124c1414c8SBarry Smith while (nz-- > 0) { /* off-diagonal elemets */ 3134c1414c8SBarry Smith /* ia[i1+1]++; */ 3144c1414c8SBarry Smith ia[i2+1]++; 3154c1414c8SBarry Smith i2++; 3164c1414c8SBarry Smith while (((col = *j++ + ishift) < tns[i2]) && nz > 0) {nz--;} 3174c1414c8SBarry Smith if (nz > 0) i2 = tvc[col]; 3184c1414c8SBarry Smith } 3194c1414c8SBarry Smith } 3204c1414c8SBarry Smith 3214c1414c8SBarry Smith /* shift ia[i] to point to next col */ 3224c1414c8SBarry Smith for (i1=1; i1<nslim_col+1; i1++) { 3234c1414c8SBarry Smith col = ia[i1-1]; 3244c1414c8SBarry Smith ia[i1] += col; 3254c1414c8SBarry Smith work[i1-1] = col - oshift; 3264c1414c8SBarry Smith } 3274c1414c8SBarry Smith 3284c1414c8SBarry Smith /* allocate space for column pointers */ 3294c1414c8SBarry Smith nz = ia[nslim_col] + (!ishift); 3304c1414c8SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&ja);CHKERRQ(ierr); 3314c1414c8SBarry Smith *jja = ja; 3324c1414c8SBarry Smith 3334c1414c8SBarry Smith /* loop over matrix putting into ja */ 3344c1414c8SBarry Smith for (i1=0,row=0; i1<nslim_row; row+=ns_row[i1],i1++) { 3354c1414c8SBarry Smith j = aj + ai[row] + ishift; 3364c1414c8SBarry Smith i2 = 0; /* Col inode index */ 3374c1414c8SBarry Smith col = *j++ + ishift; 3384c1414c8SBarry Smith i2 = tvc[col]; 3394c1414c8SBarry Smith nz = ai[row+1] - ai[row]; 3404c1414c8SBarry Smith while (nz-- > 0) { 3414c1414c8SBarry Smith /* ja[work[i1]++] = i2 + oshift; */ 3424c1414c8SBarry Smith ja[work[i2]++] = i1 + oshift; 3434c1414c8SBarry Smith i2++; 3444c1414c8SBarry Smith while(((col = *j++ + ishift) < tns[i2]) && nz > 0) {nz--;} 3454c1414c8SBarry Smith if (nz > 0) i2 = tvc[col]; 3464c1414c8SBarry Smith } 3474c1414c8SBarry Smith } 3484c1414c8SBarry Smith ierr = PetscFree(ns_col);CHKERRQ(ierr); 3494c1414c8SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 3504c1414c8SBarry Smith ierr = PetscFree(tns);CHKERRQ(ierr); 3514c1414c8SBarry Smith ierr = PetscFree(tvc);CHKERRQ(ierr); 3524c1414c8SBarry Smith PetscFunctionReturn(0); 3534c1414c8SBarry Smith } 3544c1414c8SBarry Smith 3554c1414c8SBarry Smith #undef __FUNCT__ 3564c1414c8SBarry Smith #define __FUNCT__ "MatGetColumnIJ_Inode" 3578f7157efSSatish Balay static PetscErrorCode MatGetColumnIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 3584c1414c8SBarry Smith { 3594c1414c8SBarry Smith PetscErrorCode ierr; 3604c1414c8SBarry Smith 3614c1414c8SBarry Smith PetscFunctionBegin; 3624c1414c8SBarry Smith ierr = Mat_CreateColInode(A,n,PETSC_NULL);CHKERRQ(ierr); 3634c1414c8SBarry Smith if (!ia) PetscFunctionReturn(0); 3644c1414c8SBarry Smith 3658f7157efSSatish Balay if (!blockcompressed) { 3668f7157efSSatish Balay ierr = MatGetColumnIJ_SeqAIJ(A,oshift,symmetric,blockcompressed,n,ia,ja,done);CHKERRQ(ierr);; 3678f7157efSSatish Balay } else if (symmetric) { 3684c1414c8SBarry Smith /* Since the indices are symmetric it does'nt matter */ 3694c1414c8SBarry Smith ierr = MatGetRowIJ_Inode_Symmetric(A,ia,ja,0,oshift);CHKERRQ(ierr); 3704c1414c8SBarry Smith } else { 3714c1414c8SBarry Smith ierr = MatGetColumnIJ_Inode_Nonsymmetric(A,ia,ja,0,oshift);CHKERRQ(ierr); 3724c1414c8SBarry Smith } 3734c1414c8SBarry Smith PetscFunctionReturn(0); 3744c1414c8SBarry Smith } 3754c1414c8SBarry Smith 3764c1414c8SBarry Smith #undef __FUNCT__ 3774c1414c8SBarry Smith #define __FUNCT__ "MatRestoreColumnIJ_Inode" 3788f7157efSSatish Balay static PetscErrorCode MatRestoreColumnIJ_Inode(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 3794c1414c8SBarry Smith { 3804c1414c8SBarry Smith PetscErrorCode ierr; 3814c1414c8SBarry Smith 3824c1414c8SBarry Smith PetscFunctionBegin; 3834c1414c8SBarry Smith if (!ia) PetscFunctionReturn(0); 3848f7157efSSatish Balay if (!blockcompressed) { 3858f7157efSSatish Balay ierr = MatRestoreColumnIJ_SeqAIJ(A,oshift,symmetric,blockcompressed,n,ia,ja,done);CHKERRQ(ierr);; 3868f7157efSSatish Balay } else { 3874c1414c8SBarry Smith ierr = PetscFree(*ia);CHKERRQ(ierr); 3884c1414c8SBarry Smith ierr = PetscFree(*ja);CHKERRQ(ierr); 3898f7157efSSatish Balay } 3904c1414c8SBarry Smith PetscFunctionReturn(0); 3914c1414c8SBarry Smith } 3924c1414c8SBarry Smith 3934c1414c8SBarry Smith /* ----------------------------------------------------------- */ 3944c1414c8SBarry Smith 3954c1414c8SBarry Smith #undef __FUNCT__ 3964c1414c8SBarry Smith #define __FUNCT__ "MatMult_Inode" 3974c1414c8SBarry Smith static PetscErrorCode MatMult_Inode(Mat A,Vec xx,Vec yy) 3984c1414c8SBarry Smith { 3994c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 4004c1414c8SBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,tmp0,tmp1; 401d9fead3dSBarry Smith PetscScalar *y; 402dd6ea824SBarry Smith const PetscScalar *x; 403dd6ea824SBarry Smith const MatScalar *v1,*v2,*v3,*v4,*v5; 4044c1414c8SBarry Smith PetscErrorCode ierr; 40598c9bda7SSatish Balay PetscInt *idx,i1,i2,n,i,row,node_max,*ns,*ii,nsz,sz,nonzerorow=0; 4064c1414c8SBarry Smith 4074c1414c8SBarry Smith #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 4084c1414c8SBarry Smith #pragma disjoint(*x,*y,*v1,*v2,*v3,*v4,*v5) 4094c1414c8SBarry Smith #endif 4104c1414c8SBarry Smith 4114c1414c8SBarry Smith PetscFunctionBegin; 4124c1414c8SBarry Smith if (!a->inode.size) SETERRQ(PETSC_ERR_COR,"Missing Inode Structure"); 4134c1414c8SBarry Smith node_max = a->inode.node_count; 4144c1414c8SBarry Smith ns = a->inode.size; /* Node Size array */ 415d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 4164c1414c8SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 4174c1414c8SBarry Smith idx = a->j; 4184c1414c8SBarry Smith v1 = a->a; 4194c1414c8SBarry Smith ii = a->i; 4204c1414c8SBarry Smith 4214c1414c8SBarry Smith for (i = 0,row = 0; i< node_max; ++i){ 4224c1414c8SBarry Smith nsz = ns[i]; 4234c1414c8SBarry Smith n = ii[1] - ii[0]; 42498c9bda7SSatish Balay nonzerorow += (n>0)*nsz; 4254c1414c8SBarry Smith ii += nsz; 4264c1414c8SBarry Smith sz = n; /* No of non zeros in this row */ 4274c1414c8SBarry Smith /* Switch on the size of Node */ 4284c1414c8SBarry Smith switch (nsz){ /* Each loop in 'case' is unrolled */ 4294c1414c8SBarry Smith case 1 : 4304c1414c8SBarry Smith sum1 = 0; 4314c1414c8SBarry Smith 4324c1414c8SBarry Smith for(n = 0; n< sz-1; n+=2) { 4334c1414c8SBarry Smith i1 = idx[0]; /* The instructions are ordered to */ 4344c1414c8SBarry Smith i2 = idx[1]; /* make the compiler's job easy */ 4354c1414c8SBarry Smith idx += 2; 4364c1414c8SBarry Smith tmp0 = x[i1]; 4374c1414c8SBarry Smith tmp1 = x[i2]; 4384c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 4394c1414c8SBarry Smith } 4404c1414c8SBarry Smith 4414c1414c8SBarry Smith if (n == sz-1){ /* Take care of the last nonzero */ 4424c1414c8SBarry Smith tmp0 = x[*idx++]; 4434c1414c8SBarry Smith sum1 += *v1++ * tmp0; 4444c1414c8SBarry Smith } 4454c1414c8SBarry Smith y[row++]=sum1; 4464c1414c8SBarry Smith break; 4474c1414c8SBarry Smith case 2: 4484c1414c8SBarry Smith sum1 = 0; 4494c1414c8SBarry Smith sum2 = 0; 4504c1414c8SBarry Smith v2 = v1 + n; 4514c1414c8SBarry Smith 4524c1414c8SBarry Smith for (n = 0; n< sz-1; n+=2) { 4534c1414c8SBarry Smith i1 = idx[0]; 4544c1414c8SBarry Smith i2 = idx[1]; 4554c1414c8SBarry Smith idx += 2; 4564c1414c8SBarry Smith tmp0 = x[i1]; 4574c1414c8SBarry Smith tmp1 = x[i2]; 4584c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 4594c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 4604c1414c8SBarry Smith } 4614c1414c8SBarry Smith if (n == sz-1){ 4624c1414c8SBarry Smith tmp0 = x[*idx++]; 4634c1414c8SBarry Smith sum1 += *v1++ * tmp0; 4644c1414c8SBarry Smith sum2 += *v2++ * tmp0; 4654c1414c8SBarry Smith } 4664c1414c8SBarry Smith y[row++]=sum1; 4674c1414c8SBarry Smith y[row++]=sum2; 4684c1414c8SBarry Smith v1 =v2; /* Since the next block to be processed starts there*/ 4694c1414c8SBarry Smith idx +=sz; 4704c1414c8SBarry Smith break; 4714c1414c8SBarry Smith case 3: 4724c1414c8SBarry Smith sum1 = 0; 4734c1414c8SBarry Smith sum2 = 0; 4744c1414c8SBarry Smith sum3 = 0; 4754c1414c8SBarry Smith v2 = v1 + n; 4764c1414c8SBarry Smith v3 = v2 + n; 4774c1414c8SBarry Smith 4784c1414c8SBarry Smith for (n = 0; n< sz-1; n+=2) { 4794c1414c8SBarry Smith i1 = idx[0]; 4804c1414c8SBarry Smith i2 = idx[1]; 4814c1414c8SBarry Smith idx += 2; 4824c1414c8SBarry Smith tmp0 = x[i1]; 4834c1414c8SBarry Smith tmp1 = x[i2]; 4844c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 4854c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 4864c1414c8SBarry Smith sum3 += v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 4874c1414c8SBarry Smith } 4884c1414c8SBarry Smith if (n == sz-1){ 4894c1414c8SBarry Smith tmp0 = x[*idx++]; 4904c1414c8SBarry Smith sum1 += *v1++ * tmp0; 4914c1414c8SBarry Smith sum2 += *v2++ * tmp0; 4924c1414c8SBarry Smith sum3 += *v3++ * tmp0; 4934c1414c8SBarry Smith } 4944c1414c8SBarry Smith y[row++]=sum1; 4954c1414c8SBarry Smith y[row++]=sum2; 4964c1414c8SBarry Smith y[row++]=sum3; 4974c1414c8SBarry Smith v1 =v3; /* Since the next block to be processed starts there*/ 4984c1414c8SBarry Smith idx +=2*sz; 4994c1414c8SBarry Smith break; 5004c1414c8SBarry Smith case 4: 5014c1414c8SBarry Smith sum1 = 0; 5024c1414c8SBarry Smith sum2 = 0; 5034c1414c8SBarry Smith sum3 = 0; 5044c1414c8SBarry Smith sum4 = 0; 5054c1414c8SBarry Smith v2 = v1 + n; 5064c1414c8SBarry Smith v3 = v2 + n; 5074c1414c8SBarry Smith v4 = v3 + n; 5084c1414c8SBarry Smith 5094c1414c8SBarry Smith for (n = 0; n< sz-1; n+=2) { 5104c1414c8SBarry Smith i1 = idx[0]; 5114c1414c8SBarry Smith i2 = idx[1]; 5124c1414c8SBarry Smith idx += 2; 5134c1414c8SBarry Smith tmp0 = x[i1]; 5144c1414c8SBarry Smith tmp1 = x[i2]; 5154c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2; 5164c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2; 5174c1414c8SBarry Smith sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2; 5184c1414c8SBarry Smith sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2; 5194c1414c8SBarry Smith } 5204c1414c8SBarry Smith if (n == sz-1){ 5214c1414c8SBarry Smith tmp0 = x[*idx++]; 5224c1414c8SBarry Smith sum1 += *v1++ * tmp0; 5234c1414c8SBarry Smith sum2 += *v2++ * tmp0; 5244c1414c8SBarry Smith sum3 += *v3++ * tmp0; 5254c1414c8SBarry Smith sum4 += *v4++ * tmp0; 5264c1414c8SBarry Smith } 5274c1414c8SBarry Smith y[row++]=sum1; 5284c1414c8SBarry Smith y[row++]=sum2; 5294c1414c8SBarry Smith y[row++]=sum3; 5304c1414c8SBarry Smith y[row++]=sum4; 5314c1414c8SBarry Smith v1 =v4; /* Since the next block to be processed starts there*/ 5324c1414c8SBarry Smith idx +=3*sz; 5334c1414c8SBarry Smith break; 5344c1414c8SBarry Smith case 5: 5354c1414c8SBarry Smith sum1 = 0; 5364c1414c8SBarry Smith sum2 = 0; 5374c1414c8SBarry Smith sum3 = 0; 5384c1414c8SBarry Smith sum4 = 0; 5394c1414c8SBarry Smith sum5 = 0; 5404c1414c8SBarry Smith v2 = v1 + n; 5414c1414c8SBarry Smith v3 = v2 + n; 5424c1414c8SBarry Smith v4 = v3 + n; 5434c1414c8SBarry Smith v5 = v4 + n; 5444c1414c8SBarry Smith 5454c1414c8SBarry Smith for (n = 0; n<sz-1; n+=2) { 5464c1414c8SBarry Smith i1 = idx[0]; 5474c1414c8SBarry Smith i2 = idx[1]; 5484c1414c8SBarry Smith idx += 2; 5494c1414c8SBarry Smith tmp0 = x[i1]; 5504c1414c8SBarry Smith tmp1 = x[i2]; 5514c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2; 5524c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2; 5534c1414c8SBarry Smith sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2; 5544c1414c8SBarry Smith sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2; 5554c1414c8SBarry Smith sum5 += v5[0] * tmp0 + v5[1] *tmp1; v5 += 2; 5564c1414c8SBarry Smith } 5574c1414c8SBarry Smith if (n == sz-1){ 5584c1414c8SBarry Smith tmp0 = x[*idx++]; 5594c1414c8SBarry Smith sum1 += *v1++ * tmp0; 5604c1414c8SBarry Smith sum2 += *v2++ * tmp0; 5614c1414c8SBarry Smith sum3 += *v3++ * tmp0; 5624c1414c8SBarry Smith sum4 += *v4++ * tmp0; 5634c1414c8SBarry Smith sum5 += *v5++ * tmp0; 5644c1414c8SBarry Smith } 5654c1414c8SBarry Smith y[row++]=sum1; 5664c1414c8SBarry Smith y[row++]=sum2; 5674c1414c8SBarry Smith y[row++]=sum3; 5684c1414c8SBarry Smith y[row++]=sum4; 5694c1414c8SBarry Smith y[row++]=sum5; 5704c1414c8SBarry Smith v1 =v5; /* Since the next block to be processed starts there */ 5714c1414c8SBarry Smith idx +=4*sz; 5724c1414c8SBarry Smith break; 5734c1414c8SBarry Smith default : 5744c1414c8SBarry Smith SETERRQ(PETSC_ERR_COR,"Node size not yet supported"); 5754c1414c8SBarry Smith } 5764c1414c8SBarry Smith } 577d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 5784c1414c8SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 57998c9bda7SSatish Balay ierr = PetscLogFlops(2*a->nz - nonzerorow);CHKERRQ(ierr); 5804c1414c8SBarry Smith PetscFunctionReturn(0); 5814c1414c8SBarry Smith } 5824c1414c8SBarry Smith /* ----------------------------------------------------------- */ 5834c1414c8SBarry Smith /* Almost same code as the MatMult_Inode() */ 5844c1414c8SBarry Smith #undef __FUNCT__ 5854c1414c8SBarry Smith #define __FUNCT__ "MatMultAdd_Inode" 5864c1414c8SBarry Smith static PetscErrorCode MatMultAdd_Inode(Mat A,Vec xx,Vec zz,Vec yy) 5874c1414c8SBarry Smith { 5884c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 5894c1414c8SBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,tmp0,tmp1; 590dd6ea824SBarry Smith MatScalar *v1,*v2,*v3,*v4,*v5; 591dd6ea824SBarry Smith PetscScalar *x,*y,*z,*zt; 5924c1414c8SBarry Smith PetscErrorCode ierr; 5934c1414c8SBarry Smith PetscInt *idx,i1,i2,n,i,row,node_max,*ns,*ii,nsz,sz; 5944c1414c8SBarry Smith 5954c1414c8SBarry Smith PetscFunctionBegin; 5964c1414c8SBarry Smith if (!a->inode.size) SETERRQ(PETSC_ERR_COR,"Missing Inode Structure"); 5974c1414c8SBarry Smith node_max = a->inode.node_count; 5984c1414c8SBarry Smith ns = a->inode.size; /* Node Size array */ 5994c1414c8SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6004c1414c8SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6014c1414c8SBarry Smith if (zz != yy) { 6024c1414c8SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 6034c1414c8SBarry Smith } else { 6044c1414c8SBarry Smith z = y; 6054c1414c8SBarry Smith } 6064c1414c8SBarry Smith zt = z; 6074c1414c8SBarry Smith 6084c1414c8SBarry Smith idx = a->j; 6094c1414c8SBarry Smith v1 = a->a; 6104c1414c8SBarry Smith ii = a->i; 6114c1414c8SBarry Smith 6124c1414c8SBarry Smith for (i = 0,row = 0; i< node_max; ++i){ 6134c1414c8SBarry Smith nsz = ns[i]; 6144c1414c8SBarry Smith n = ii[1] - ii[0]; 6154c1414c8SBarry Smith ii += nsz; 6164c1414c8SBarry Smith sz = n; /* No of non zeros in this row */ 6174c1414c8SBarry Smith /* Switch on the size of Node */ 6184c1414c8SBarry Smith switch (nsz){ /* Each loop in 'case' is unrolled */ 6194c1414c8SBarry Smith case 1 : 6204c1414c8SBarry Smith sum1 = *zt++; 6214c1414c8SBarry Smith 6224c1414c8SBarry Smith for(n = 0; n< sz-1; n+=2) { 6234c1414c8SBarry Smith i1 = idx[0]; /* The instructions are ordered to */ 6244c1414c8SBarry Smith i2 = idx[1]; /* make the compiler's job easy */ 6254c1414c8SBarry Smith idx += 2; 6264c1414c8SBarry Smith tmp0 = x[i1]; 6274c1414c8SBarry Smith tmp1 = x[i2]; 6284c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 6294c1414c8SBarry Smith } 6304c1414c8SBarry Smith 6314c1414c8SBarry Smith if(n == sz-1){ /* Take care of the last nonzero */ 6324c1414c8SBarry Smith tmp0 = x[*idx++]; 6334c1414c8SBarry Smith sum1 += *v1++ * tmp0; 6344c1414c8SBarry Smith } 6354c1414c8SBarry Smith y[row++]=sum1; 6364c1414c8SBarry Smith break; 6374c1414c8SBarry Smith case 2: 6384c1414c8SBarry Smith sum1 = *zt++; 6394c1414c8SBarry Smith sum2 = *zt++; 6404c1414c8SBarry Smith v2 = v1 + n; 6414c1414c8SBarry Smith 6424c1414c8SBarry Smith for(n = 0; n< sz-1; n+=2) { 6434c1414c8SBarry Smith i1 = idx[0]; 6444c1414c8SBarry Smith i2 = idx[1]; 6454c1414c8SBarry Smith idx += 2; 6464c1414c8SBarry Smith tmp0 = x[i1]; 6474c1414c8SBarry Smith tmp1 = x[i2]; 6484c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 6494c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 6504c1414c8SBarry Smith } 6514c1414c8SBarry Smith if(n == sz-1){ 6524c1414c8SBarry Smith tmp0 = x[*idx++]; 6534c1414c8SBarry Smith sum1 += *v1++ * tmp0; 6544c1414c8SBarry Smith sum2 += *v2++ * tmp0; 6554c1414c8SBarry Smith } 6564c1414c8SBarry Smith y[row++]=sum1; 6574c1414c8SBarry Smith y[row++]=sum2; 6584c1414c8SBarry Smith v1 =v2; /* Since the next block to be processed starts there*/ 6594c1414c8SBarry Smith idx +=sz; 6604c1414c8SBarry Smith break; 6614c1414c8SBarry Smith case 3: 6624c1414c8SBarry Smith sum1 = *zt++; 6634c1414c8SBarry Smith sum2 = *zt++; 6644c1414c8SBarry Smith sum3 = *zt++; 6654c1414c8SBarry Smith v2 = v1 + n; 6664c1414c8SBarry Smith v3 = v2 + n; 6674c1414c8SBarry Smith 6684c1414c8SBarry Smith for (n = 0; n< sz-1; n+=2) { 6694c1414c8SBarry Smith i1 = idx[0]; 6704c1414c8SBarry Smith i2 = idx[1]; 6714c1414c8SBarry Smith idx += 2; 6724c1414c8SBarry Smith tmp0 = x[i1]; 6734c1414c8SBarry Smith tmp1 = x[i2]; 6744c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 6754c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 6764c1414c8SBarry Smith sum3 += v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 6774c1414c8SBarry Smith } 6784c1414c8SBarry Smith if (n == sz-1){ 6794c1414c8SBarry Smith tmp0 = x[*idx++]; 6804c1414c8SBarry Smith sum1 += *v1++ * tmp0; 6814c1414c8SBarry Smith sum2 += *v2++ * tmp0; 6824c1414c8SBarry Smith sum3 += *v3++ * tmp0; 6834c1414c8SBarry Smith } 6844c1414c8SBarry Smith y[row++]=sum1; 6854c1414c8SBarry Smith y[row++]=sum2; 6864c1414c8SBarry Smith y[row++]=sum3; 6874c1414c8SBarry Smith v1 =v3; /* Since the next block to be processed starts there*/ 6884c1414c8SBarry Smith idx +=2*sz; 6894c1414c8SBarry Smith break; 6904c1414c8SBarry Smith case 4: 6914c1414c8SBarry Smith sum1 = *zt++; 6924c1414c8SBarry Smith sum2 = *zt++; 6934c1414c8SBarry Smith sum3 = *zt++; 6944c1414c8SBarry Smith sum4 = *zt++; 6954c1414c8SBarry Smith v2 = v1 + n; 6964c1414c8SBarry Smith v3 = v2 + n; 6974c1414c8SBarry Smith v4 = v3 + n; 6984c1414c8SBarry Smith 6994c1414c8SBarry Smith for (n = 0; n< sz-1; n+=2) { 7004c1414c8SBarry Smith i1 = idx[0]; 7014c1414c8SBarry Smith i2 = idx[1]; 7024c1414c8SBarry Smith idx += 2; 7034c1414c8SBarry Smith tmp0 = x[i1]; 7044c1414c8SBarry Smith tmp1 = x[i2]; 7054c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2; 7064c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2; 7074c1414c8SBarry Smith sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2; 7084c1414c8SBarry Smith sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2; 7094c1414c8SBarry Smith } 7104c1414c8SBarry Smith if (n == sz-1){ 7114c1414c8SBarry Smith tmp0 = x[*idx++]; 7124c1414c8SBarry Smith sum1 += *v1++ * tmp0; 7134c1414c8SBarry Smith sum2 += *v2++ * tmp0; 7144c1414c8SBarry Smith sum3 += *v3++ * tmp0; 7154c1414c8SBarry Smith sum4 += *v4++ * tmp0; 7164c1414c8SBarry Smith } 7174c1414c8SBarry Smith y[row++]=sum1; 7184c1414c8SBarry Smith y[row++]=sum2; 7194c1414c8SBarry Smith y[row++]=sum3; 7204c1414c8SBarry Smith y[row++]=sum4; 7214c1414c8SBarry Smith v1 =v4; /* Since the next block to be processed starts there*/ 7224c1414c8SBarry Smith idx +=3*sz; 7234c1414c8SBarry Smith break; 7244c1414c8SBarry Smith case 5: 7254c1414c8SBarry Smith sum1 = *zt++; 7264c1414c8SBarry Smith sum2 = *zt++; 7274c1414c8SBarry Smith sum3 = *zt++; 7284c1414c8SBarry Smith sum4 = *zt++; 7294c1414c8SBarry Smith sum5 = *zt++; 7304c1414c8SBarry Smith v2 = v1 + n; 7314c1414c8SBarry Smith v3 = v2 + n; 7324c1414c8SBarry Smith v4 = v3 + n; 7334c1414c8SBarry Smith v5 = v4 + n; 7344c1414c8SBarry Smith 7354c1414c8SBarry Smith for (n = 0; n<sz-1; n+=2) { 7364c1414c8SBarry Smith i1 = idx[0]; 7374c1414c8SBarry Smith i2 = idx[1]; 7384c1414c8SBarry Smith idx += 2; 7394c1414c8SBarry Smith tmp0 = x[i1]; 7404c1414c8SBarry Smith tmp1 = x[i2]; 7414c1414c8SBarry Smith sum1 += v1[0] * tmp0 + v1[1] *tmp1; v1 += 2; 7424c1414c8SBarry Smith sum2 += v2[0] * tmp0 + v2[1] *tmp1; v2 += 2; 7434c1414c8SBarry Smith sum3 += v3[0] * tmp0 + v3[1] *tmp1; v3 += 2; 7444c1414c8SBarry Smith sum4 += v4[0] * tmp0 + v4[1] *tmp1; v4 += 2; 7454c1414c8SBarry Smith sum5 += v5[0] * tmp0 + v5[1] *tmp1; v5 += 2; 7464c1414c8SBarry Smith } 7474c1414c8SBarry Smith if(n == sz-1){ 7484c1414c8SBarry Smith tmp0 = x[*idx++]; 7494c1414c8SBarry Smith sum1 += *v1++ * tmp0; 7504c1414c8SBarry Smith sum2 += *v2++ * tmp0; 7514c1414c8SBarry Smith sum3 += *v3++ * tmp0; 7524c1414c8SBarry Smith sum4 += *v4++ * tmp0; 7534c1414c8SBarry Smith sum5 += *v5++ * tmp0; 7544c1414c8SBarry Smith } 7554c1414c8SBarry Smith y[row++]=sum1; 7564c1414c8SBarry Smith y[row++]=sum2; 7574c1414c8SBarry Smith y[row++]=sum3; 7584c1414c8SBarry Smith y[row++]=sum4; 7594c1414c8SBarry Smith y[row++]=sum5; 7604c1414c8SBarry Smith v1 =v5; /* Since the next block to be processed starts there */ 7614c1414c8SBarry Smith idx +=4*sz; 7624c1414c8SBarry Smith break; 7634c1414c8SBarry Smith default : 7644c1414c8SBarry Smith SETERRQ(PETSC_ERR_COR,"Node size not yet supported"); 7654c1414c8SBarry Smith } 7664c1414c8SBarry Smith } 7674c1414c8SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 7684c1414c8SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 7694c1414c8SBarry Smith if (zz != yy) { 7704c1414c8SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 7714c1414c8SBarry Smith } 7724c1414c8SBarry Smith ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 7734c1414c8SBarry Smith PetscFunctionReturn(0); 7744c1414c8SBarry Smith } 7754c1414c8SBarry Smith 7764c1414c8SBarry Smith /* ----------------------------------------------------------- */ 7774c1414c8SBarry Smith #undef __FUNCT__ 7784c1414c8SBarry Smith #define __FUNCT__ "MatSolve_Inode" 7794c1414c8SBarry Smith PetscErrorCode MatSolve_Inode(Mat A,Vec bb,Vec xx) 7804c1414c8SBarry Smith { 7814c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7824c1414c8SBarry Smith IS iscol = a->col,isrow = a->row; 7834c1414c8SBarry Smith PetscErrorCode ierr; 784d0f46423SBarry Smith PetscInt *r,*c,i,j,n = A->rmap->n,*ai = a->i,nz,*a_j = a->j; 7854c1414c8SBarry Smith PetscInt node_max,*ns,row,nsz,aii,*vi,*ad,*aj,i0,i1,*rout,*cout; 786d9fead3dSBarry Smith PetscScalar *x,*tmp,*tmps,tmp0,tmp1; 787d9fead3dSBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5; 788dd6ea824SBarry Smith const MatScalar *v1,*v2,*v3,*v4,*v5,*a_a = a->a,*aa; 789dd6ea824SBarry Smith const PetscScalar *b; 7904c1414c8SBarry Smith 7914c1414c8SBarry Smith PetscFunctionBegin; 7924c1414c8SBarry Smith if (!a->inode.size) SETERRQ(PETSC_ERR_COR,"Missing Inode Structure"); 7934c1414c8SBarry Smith node_max = a->inode.node_count; 7944c1414c8SBarry Smith ns = a->inode.size; /* Node Size array */ 7954c1414c8SBarry Smith 796d9fead3dSBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 7974c1414c8SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 7984c1414c8SBarry Smith tmp = a->solve_work; 7994c1414c8SBarry Smith 8004c1414c8SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8014c1414c8SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 8024c1414c8SBarry Smith 8034c1414c8SBarry Smith /* forward solve the lower triangular */ 8044c1414c8SBarry Smith tmps = tmp ; 8054c1414c8SBarry Smith aa = a_a ; 8064c1414c8SBarry Smith aj = a_j ; 8074c1414c8SBarry Smith ad = a->diag; 8084c1414c8SBarry Smith 8094c1414c8SBarry Smith for (i = 0,row = 0; i< node_max; ++i){ 8104c1414c8SBarry Smith nsz = ns[i]; 8114c1414c8SBarry Smith aii = ai[row]; 8124c1414c8SBarry Smith v1 = aa + aii; 8134c1414c8SBarry Smith vi = aj + aii; 8144c1414c8SBarry Smith nz = ad[row]- aii; 8154c1414c8SBarry Smith 8164c1414c8SBarry Smith switch (nsz){ /* Each loop in 'case' is unrolled */ 8174c1414c8SBarry Smith case 1 : 8184c1414c8SBarry Smith sum1 = b[*r++]; 8194c1414c8SBarry Smith /* while (nz--) sum1 -= *v1++ *tmps[*vi++];*/ 8204c1414c8SBarry Smith for(j=0; j<nz-1; j+=2){ 8214c1414c8SBarry Smith i0 = vi[0]; 8224c1414c8SBarry Smith i1 = vi[1]; 8234c1414c8SBarry Smith vi +=2; 8244c1414c8SBarry Smith tmp0 = tmps[i0]; 8254c1414c8SBarry Smith tmp1 = tmps[i1]; 8264c1414c8SBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 8274c1414c8SBarry Smith } 8284c1414c8SBarry Smith if(j == nz-1){ 8294c1414c8SBarry Smith tmp0 = tmps[*vi++]; 8304c1414c8SBarry Smith sum1 -= *v1++ *tmp0; 8314c1414c8SBarry Smith } 8324c1414c8SBarry Smith tmp[row ++]=sum1; 8334c1414c8SBarry Smith break; 8344c1414c8SBarry Smith case 2: 8354c1414c8SBarry Smith sum1 = b[*r++]; 8364c1414c8SBarry Smith sum2 = b[*r++]; 8374c1414c8SBarry Smith v2 = aa + ai[row+1]; 8384c1414c8SBarry Smith 8394c1414c8SBarry Smith for(j=0; j<nz-1; j+=2){ 8404c1414c8SBarry Smith i0 = vi[0]; 8414c1414c8SBarry Smith i1 = vi[1]; 8424c1414c8SBarry Smith vi +=2; 8434c1414c8SBarry Smith tmp0 = tmps[i0]; 8444c1414c8SBarry Smith tmp1 = tmps[i1]; 8454c1414c8SBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 8464c1414c8SBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 8474c1414c8SBarry Smith } 8484c1414c8SBarry Smith if(j == nz-1){ 8494c1414c8SBarry Smith tmp0 = tmps[*vi++]; 8504c1414c8SBarry Smith sum1 -= *v1++ *tmp0; 8514c1414c8SBarry Smith sum2 -= *v2++ *tmp0; 8524c1414c8SBarry Smith } 8534c1414c8SBarry Smith sum2 -= *v2++ * sum1; 8544c1414c8SBarry Smith tmp[row ++]=sum1; 8554c1414c8SBarry Smith tmp[row ++]=sum2; 8564c1414c8SBarry Smith break; 8574c1414c8SBarry Smith case 3: 8584c1414c8SBarry Smith sum1 = b[*r++]; 8594c1414c8SBarry Smith sum2 = b[*r++]; 8604c1414c8SBarry Smith sum3 = b[*r++]; 8614c1414c8SBarry Smith v2 = aa + ai[row+1]; 8624c1414c8SBarry Smith v3 = aa + ai[row+2]; 8634c1414c8SBarry Smith 8644c1414c8SBarry Smith for (j=0; j<nz-1; j+=2){ 8654c1414c8SBarry Smith i0 = vi[0]; 8664c1414c8SBarry Smith i1 = vi[1]; 8674c1414c8SBarry Smith vi +=2; 8684c1414c8SBarry Smith tmp0 = tmps[i0]; 8694c1414c8SBarry Smith tmp1 = tmps[i1]; 8704c1414c8SBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 8714c1414c8SBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 8724c1414c8SBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 8734c1414c8SBarry Smith } 8744c1414c8SBarry Smith if (j == nz-1){ 8754c1414c8SBarry Smith tmp0 = tmps[*vi++]; 8764c1414c8SBarry Smith sum1 -= *v1++ *tmp0; 8774c1414c8SBarry Smith sum2 -= *v2++ *tmp0; 8784c1414c8SBarry Smith sum3 -= *v3++ *tmp0; 8794c1414c8SBarry Smith } 8804c1414c8SBarry Smith sum2 -= *v2++ * sum1; 8814c1414c8SBarry Smith sum3 -= *v3++ * sum1; 8824c1414c8SBarry Smith sum3 -= *v3++ * sum2; 8834c1414c8SBarry Smith tmp[row ++]=sum1; 8844c1414c8SBarry Smith tmp[row ++]=sum2; 8854c1414c8SBarry Smith tmp[row ++]=sum3; 8864c1414c8SBarry Smith break; 8874c1414c8SBarry Smith 8884c1414c8SBarry Smith case 4: 8894c1414c8SBarry Smith sum1 = b[*r++]; 8904c1414c8SBarry Smith sum2 = b[*r++]; 8914c1414c8SBarry Smith sum3 = b[*r++]; 8924c1414c8SBarry Smith sum4 = b[*r++]; 8934c1414c8SBarry Smith v2 = aa + ai[row+1]; 8944c1414c8SBarry Smith v3 = aa + ai[row+2]; 8954c1414c8SBarry Smith v4 = aa + ai[row+3]; 8964c1414c8SBarry Smith 8974c1414c8SBarry Smith for (j=0; j<nz-1; j+=2){ 8984c1414c8SBarry Smith i0 = vi[0]; 8994c1414c8SBarry Smith i1 = vi[1]; 9004c1414c8SBarry Smith vi +=2; 9014c1414c8SBarry Smith tmp0 = tmps[i0]; 9024c1414c8SBarry Smith tmp1 = tmps[i1]; 9034c1414c8SBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 9044c1414c8SBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 9054c1414c8SBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 9064c1414c8SBarry Smith sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2; 9074c1414c8SBarry Smith } 9084c1414c8SBarry Smith if (j == nz-1){ 9094c1414c8SBarry Smith tmp0 = tmps[*vi++]; 9104c1414c8SBarry Smith sum1 -= *v1++ *tmp0; 9114c1414c8SBarry Smith sum2 -= *v2++ *tmp0; 9124c1414c8SBarry Smith sum3 -= *v3++ *tmp0; 9134c1414c8SBarry Smith sum4 -= *v4++ *tmp0; 9144c1414c8SBarry Smith } 9154c1414c8SBarry Smith sum2 -= *v2++ * sum1; 9164c1414c8SBarry Smith sum3 -= *v3++ * sum1; 9174c1414c8SBarry Smith sum4 -= *v4++ * sum1; 9184c1414c8SBarry Smith sum3 -= *v3++ * sum2; 9194c1414c8SBarry Smith sum4 -= *v4++ * sum2; 9204c1414c8SBarry Smith sum4 -= *v4++ * sum3; 9214c1414c8SBarry Smith 9224c1414c8SBarry Smith tmp[row ++]=sum1; 9234c1414c8SBarry Smith tmp[row ++]=sum2; 9244c1414c8SBarry Smith tmp[row ++]=sum3; 9254c1414c8SBarry Smith tmp[row ++]=sum4; 9264c1414c8SBarry Smith break; 9274c1414c8SBarry Smith case 5: 9284c1414c8SBarry Smith sum1 = b[*r++]; 9294c1414c8SBarry Smith sum2 = b[*r++]; 9304c1414c8SBarry Smith sum3 = b[*r++]; 9314c1414c8SBarry Smith sum4 = b[*r++]; 9324c1414c8SBarry Smith sum5 = b[*r++]; 9334c1414c8SBarry Smith v2 = aa + ai[row+1]; 9344c1414c8SBarry Smith v3 = aa + ai[row+2]; 9354c1414c8SBarry Smith v4 = aa + ai[row+3]; 9364c1414c8SBarry Smith v5 = aa + ai[row+4]; 9374c1414c8SBarry Smith 9384c1414c8SBarry Smith for (j=0; j<nz-1; j+=2){ 9394c1414c8SBarry Smith i0 = vi[0]; 9404c1414c8SBarry Smith i1 = vi[1]; 9414c1414c8SBarry Smith vi +=2; 9424c1414c8SBarry Smith tmp0 = tmps[i0]; 9434c1414c8SBarry Smith tmp1 = tmps[i1]; 9444c1414c8SBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 9454c1414c8SBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 9464c1414c8SBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 9474c1414c8SBarry Smith sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2; 9484c1414c8SBarry Smith sum5 -= v5[0] * tmp0 + v5[1] * tmp1; v5 += 2; 9494c1414c8SBarry Smith } 9504c1414c8SBarry Smith if (j == nz-1){ 9514c1414c8SBarry Smith tmp0 = tmps[*vi++]; 9524c1414c8SBarry Smith sum1 -= *v1++ *tmp0; 9534c1414c8SBarry Smith sum2 -= *v2++ *tmp0; 9544c1414c8SBarry Smith sum3 -= *v3++ *tmp0; 9554c1414c8SBarry Smith sum4 -= *v4++ *tmp0; 9564c1414c8SBarry Smith sum5 -= *v5++ *tmp0; 9574c1414c8SBarry Smith } 9584c1414c8SBarry Smith 9594c1414c8SBarry Smith sum2 -= *v2++ * sum1; 9604c1414c8SBarry Smith sum3 -= *v3++ * sum1; 9614c1414c8SBarry Smith sum4 -= *v4++ * sum1; 9624c1414c8SBarry Smith sum5 -= *v5++ * sum1; 9634c1414c8SBarry Smith sum3 -= *v3++ * sum2; 9644c1414c8SBarry Smith sum4 -= *v4++ * sum2; 9654c1414c8SBarry Smith sum5 -= *v5++ * sum2; 9664c1414c8SBarry Smith sum4 -= *v4++ * sum3; 9674c1414c8SBarry Smith sum5 -= *v5++ * sum3; 9684c1414c8SBarry Smith sum5 -= *v5++ * sum4; 9694c1414c8SBarry Smith 9704c1414c8SBarry Smith tmp[row ++]=sum1; 9714c1414c8SBarry Smith tmp[row ++]=sum2; 9724c1414c8SBarry Smith tmp[row ++]=sum3; 9734c1414c8SBarry Smith tmp[row ++]=sum4; 9744c1414c8SBarry Smith tmp[row ++]=sum5; 9754c1414c8SBarry Smith break; 9764c1414c8SBarry Smith default: 9774c1414c8SBarry Smith SETERRQ(PETSC_ERR_COR,"Node size not yet supported \n"); 9784c1414c8SBarry Smith } 9794c1414c8SBarry Smith } 9804c1414c8SBarry Smith /* backward solve the upper triangular */ 9814c1414c8SBarry Smith for (i=node_max -1 ,row = n-1 ; i>=0; i--){ 9824c1414c8SBarry Smith nsz = ns[i]; 9834c1414c8SBarry Smith aii = ai[row+1] -1; 9844c1414c8SBarry Smith v1 = aa + aii; 9854c1414c8SBarry Smith vi = aj + aii; 9864c1414c8SBarry Smith nz = aii- ad[row]; 9874c1414c8SBarry Smith switch (nsz){ /* Each loop in 'case' is unrolled */ 9884c1414c8SBarry Smith case 1 : 9894c1414c8SBarry Smith sum1 = tmp[row]; 9904c1414c8SBarry Smith 9914c1414c8SBarry Smith for(j=nz ; j>1; j-=2){ 9924c1414c8SBarry Smith vi -=2; 9934c1414c8SBarry Smith i0 = vi[2]; 9944c1414c8SBarry Smith i1 = vi[1]; 9954c1414c8SBarry Smith tmp0 = tmps[i0]; 9964c1414c8SBarry Smith tmp1 = tmps[i1]; 9974c1414c8SBarry Smith v1 -= 2; 9984c1414c8SBarry Smith sum1 -= v1[2] * tmp0 + v1[1] * tmp1; 9994c1414c8SBarry Smith } 10004c1414c8SBarry Smith if (j==1){ 10014c1414c8SBarry Smith tmp0 = tmps[*vi--]; 10024c1414c8SBarry Smith sum1 -= *v1-- * tmp0; 10034c1414c8SBarry Smith } 10044c1414c8SBarry Smith x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--; 10054c1414c8SBarry Smith break; 10064c1414c8SBarry Smith case 2 : 10074c1414c8SBarry Smith sum1 = tmp[row]; 10084c1414c8SBarry Smith sum2 = tmp[row -1]; 10094c1414c8SBarry Smith v2 = aa + ai[row]-1; 10104c1414c8SBarry Smith for (j=nz ; j>1; j-=2){ 10114c1414c8SBarry Smith vi -=2; 10124c1414c8SBarry Smith i0 = vi[2]; 10134c1414c8SBarry Smith i1 = vi[1]; 10144c1414c8SBarry Smith tmp0 = tmps[i0]; 10154c1414c8SBarry Smith tmp1 = tmps[i1]; 10164c1414c8SBarry Smith v1 -= 2; 10174c1414c8SBarry Smith v2 -= 2; 10184c1414c8SBarry Smith sum1 -= v1[2] * tmp0 + v1[1] * tmp1; 10194c1414c8SBarry Smith sum2 -= v2[2] * tmp0 + v2[1] * tmp1; 10204c1414c8SBarry Smith } 10214c1414c8SBarry Smith if (j==1){ 10224c1414c8SBarry Smith tmp0 = tmps[*vi--]; 10234c1414c8SBarry Smith sum1 -= *v1-- * tmp0; 10244c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 10254c1414c8SBarry Smith } 10264c1414c8SBarry Smith 10274c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--; 10284c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 10294c1414c8SBarry Smith x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--; 10304c1414c8SBarry Smith break; 10314c1414c8SBarry Smith case 3 : 10324c1414c8SBarry Smith sum1 = tmp[row]; 10334c1414c8SBarry Smith sum2 = tmp[row -1]; 10344c1414c8SBarry Smith sum3 = tmp[row -2]; 10354c1414c8SBarry Smith v2 = aa + ai[row]-1; 10364c1414c8SBarry Smith v3 = aa + ai[row -1]-1; 10374c1414c8SBarry Smith for (j=nz ; j>1; j-=2){ 10384c1414c8SBarry Smith vi -=2; 10394c1414c8SBarry Smith i0 = vi[2]; 10404c1414c8SBarry Smith i1 = vi[1]; 10414c1414c8SBarry Smith tmp0 = tmps[i0]; 10424c1414c8SBarry Smith tmp1 = tmps[i1]; 10434c1414c8SBarry Smith v1 -= 2; 10444c1414c8SBarry Smith v2 -= 2; 10454c1414c8SBarry Smith v3 -= 2; 10464c1414c8SBarry Smith sum1 -= v1[2] * tmp0 + v1[1] * tmp1; 10474c1414c8SBarry Smith sum2 -= v2[2] * tmp0 + v2[1] * tmp1; 10484c1414c8SBarry Smith sum3 -= v3[2] * tmp0 + v3[1] * tmp1; 10494c1414c8SBarry Smith } 10504c1414c8SBarry Smith if (j==1){ 10514c1414c8SBarry Smith tmp0 = tmps[*vi--]; 10524c1414c8SBarry Smith sum1 -= *v1-- * tmp0; 10534c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 10544c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 10554c1414c8SBarry Smith } 10564c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--; 10574c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 10584c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 10594c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--; 10604c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 10614c1414c8SBarry Smith x[*c--] = tmp[row] = sum3*a_a[ad[row]]; row--; 10624c1414c8SBarry Smith 10634c1414c8SBarry Smith break; 10644c1414c8SBarry Smith case 4 : 10654c1414c8SBarry Smith sum1 = tmp[row]; 10664c1414c8SBarry Smith sum2 = tmp[row -1]; 10674c1414c8SBarry Smith sum3 = tmp[row -2]; 10684c1414c8SBarry Smith sum4 = tmp[row -3]; 10694c1414c8SBarry Smith v2 = aa + ai[row]-1; 10704c1414c8SBarry Smith v3 = aa + ai[row -1]-1; 10714c1414c8SBarry Smith v4 = aa + ai[row -2]-1; 10724c1414c8SBarry Smith 10734c1414c8SBarry Smith for (j=nz ; j>1; j-=2){ 10744c1414c8SBarry Smith vi -=2; 10754c1414c8SBarry Smith i0 = vi[2]; 10764c1414c8SBarry Smith i1 = vi[1]; 10774c1414c8SBarry Smith tmp0 = tmps[i0]; 10784c1414c8SBarry Smith tmp1 = tmps[i1]; 10794c1414c8SBarry Smith v1 -= 2; 10804c1414c8SBarry Smith v2 -= 2; 10814c1414c8SBarry Smith v3 -= 2; 10824c1414c8SBarry Smith v4 -= 2; 10834c1414c8SBarry Smith sum1 -= v1[2] * tmp0 + v1[1] * tmp1; 10844c1414c8SBarry Smith sum2 -= v2[2] * tmp0 + v2[1] * tmp1; 10854c1414c8SBarry Smith sum3 -= v3[2] * tmp0 + v3[1] * tmp1; 10864c1414c8SBarry Smith sum4 -= v4[2] * tmp0 + v4[1] * tmp1; 10874c1414c8SBarry Smith } 10884c1414c8SBarry Smith if (j==1){ 10894c1414c8SBarry Smith tmp0 = tmps[*vi--]; 10904c1414c8SBarry Smith sum1 -= *v1-- * tmp0; 10914c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 10924c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 10934c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 10944c1414c8SBarry Smith } 10954c1414c8SBarry Smith 10964c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--; 10974c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 10984c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 10994c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11004c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--; 11014c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 11024c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11034c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum3*a_a[ad[row]]; row--; 11044c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11054c1414c8SBarry Smith x[*c--] = tmp[row] = sum4*a_a[ad[row]]; row--; 11064c1414c8SBarry Smith break; 11074c1414c8SBarry Smith case 5 : 11084c1414c8SBarry Smith sum1 = tmp[row]; 11094c1414c8SBarry Smith sum2 = tmp[row -1]; 11104c1414c8SBarry Smith sum3 = tmp[row -2]; 11114c1414c8SBarry Smith sum4 = tmp[row -3]; 11124c1414c8SBarry Smith sum5 = tmp[row -4]; 11134c1414c8SBarry Smith v2 = aa + ai[row]-1; 11144c1414c8SBarry Smith v3 = aa + ai[row -1]-1; 11154c1414c8SBarry Smith v4 = aa + ai[row -2]-1; 11164c1414c8SBarry Smith v5 = aa + ai[row -3]-1; 11174c1414c8SBarry Smith for (j=nz ; j>1; j-=2){ 11184c1414c8SBarry Smith vi -= 2; 11194c1414c8SBarry Smith i0 = vi[2]; 11204c1414c8SBarry Smith i1 = vi[1]; 11214c1414c8SBarry Smith tmp0 = tmps[i0]; 11224c1414c8SBarry Smith tmp1 = tmps[i1]; 11234c1414c8SBarry Smith v1 -= 2; 11244c1414c8SBarry Smith v2 -= 2; 11254c1414c8SBarry Smith v3 -= 2; 11264c1414c8SBarry Smith v4 -= 2; 11274c1414c8SBarry Smith v5 -= 2; 11284c1414c8SBarry Smith sum1 -= v1[2] * tmp0 + v1[1] * tmp1; 11294c1414c8SBarry Smith sum2 -= v2[2] * tmp0 + v2[1] * tmp1; 11304c1414c8SBarry Smith sum3 -= v3[2] * tmp0 + v3[1] * tmp1; 11314c1414c8SBarry Smith sum4 -= v4[2] * tmp0 + v4[1] * tmp1; 11324c1414c8SBarry Smith sum5 -= v5[2] * tmp0 + v5[1] * tmp1; 11334c1414c8SBarry Smith } 11344c1414c8SBarry Smith if (j==1){ 11354c1414c8SBarry Smith tmp0 = tmps[*vi--]; 11364c1414c8SBarry Smith sum1 -= *v1-- * tmp0; 11374c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 11384c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 11394c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11404c1414c8SBarry Smith sum5 -= *v5-- * tmp0; 11414c1414c8SBarry Smith } 11424c1414c8SBarry Smith 11434c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum1*a_a[ad[row]]; row--; 11444c1414c8SBarry Smith sum2 -= *v2-- * tmp0; 11454c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 11464c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11474c1414c8SBarry Smith sum5 -= *v5-- * tmp0; 11484c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum2*a_a[ad[row]]; row--; 11494c1414c8SBarry Smith sum3 -= *v3-- * tmp0; 11504c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11514c1414c8SBarry Smith sum5 -= *v5-- * tmp0; 11524c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum3*a_a[ad[row]]; row--; 11534c1414c8SBarry Smith sum4 -= *v4-- * tmp0; 11544c1414c8SBarry Smith sum5 -= *v5-- * tmp0; 11554c1414c8SBarry Smith tmp0 = x[*c--] = tmp[row] = sum4*a_a[ad[row]]; row--; 11564c1414c8SBarry Smith sum5 -= *v5-- * tmp0; 11574c1414c8SBarry Smith x[*c--] = tmp[row] = sum5*a_a[ad[row]]; row--; 11584c1414c8SBarry Smith break; 11594c1414c8SBarry Smith default: 11604c1414c8SBarry Smith SETERRQ(PETSC_ERR_COR,"Node size not yet supported \n"); 11614c1414c8SBarry Smith } 11624c1414c8SBarry Smith } 11634c1414c8SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 11644c1414c8SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1165d9fead3dSBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 11664c1414c8SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1167d0f46423SBarry Smith ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr); 11684c1414c8SBarry Smith PetscFunctionReturn(0); 11694c1414c8SBarry Smith } 11704c1414c8SBarry Smith 11714c1414c8SBarry Smith #undef __FUNCT__ 11724c1414c8SBarry Smith #define __FUNCT__ "MatLUFactorNumeric_Inode" 1173*719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_Inode(Mat B,Mat A,MatFactorInfo *info) 11744c1414c8SBarry Smith { 1175*719d5645SBarry Smith Mat C = B; 11764c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)C->data; 11774c1414c8SBarry Smith IS iscol = b->col,isrow = b->row,isicol = b->icol; 11784c1414c8SBarry Smith PetscErrorCode ierr; 1179d0f46423SBarry Smith PetscInt *r,*ic,*c,n = A->rmap->n,*bi = b->i; 11804c1414c8SBarry Smith PetscInt *bj = b->j,*nbj=b->j +1,*ajtmp,*bjtmp,nz,nz_tmp,row,prow; 11814c1414c8SBarry Smith PetscInt *ics,i,j,idx,*ai = a->i,*aj = a->j,*bd = b->diag,node_max,nodesz; 11824c1414c8SBarry Smith PetscInt *ns,*tmp_vec1,*tmp_vec2,*nsmap,*pj; 1183dd6ea824SBarry Smith PetscScalar mul1,mul2,mul3,tmp; 1184dd6ea824SBarry Smith MatScalar *pc1,*pc2,*pc3,*ba = b->a,*pv,*rtmp11,*rtmp22,*rtmp33; 1185dd6ea824SBarry Smith const MatScalar *v1,*v2,*v3,*aa = a->a,*rtmp1; 11864c1414c8SBarry Smith PetscReal rs=0.0; 11874c1414c8SBarry Smith LUShift_Ctx sctx; 11884c1414c8SBarry Smith PetscInt newshift; 11894c1414c8SBarry Smith 11904c1414c8SBarry Smith PetscFunctionBegin; 11914c1414c8SBarry Smith sctx.shift_top = 0; 11924c1414c8SBarry Smith sctx.nshift_max = 0; 11934c1414c8SBarry Smith sctx.shift_lo = 0; 11944c1414c8SBarry Smith sctx.shift_hi = 0; 11954c1414c8SBarry Smith 11964c1414c8SBarry Smith /* if both shift schemes are chosen by user, only use info->shiftpd */ 11974c1414c8SBarry Smith if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0; 11984c1414c8SBarry Smith if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 11994c1414c8SBarry Smith sctx.shift_top = 0; 12004c1414c8SBarry Smith for (i=0; i<n; i++) { 12014c1414c8SBarry Smith /* calculate rs = sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 12024c1414c8SBarry Smith rs = 0.0; 12034c1414c8SBarry Smith ajtmp = aj + ai[i]; 12044c1414c8SBarry Smith rtmp1 = aa + ai[i]; 12054c1414c8SBarry Smith nz = ai[i+1] - ai[i]; 12064c1414c8SBarry Smith for (j=0; j<nz; j++){ 12074c1414c8SBarry Smith if (*ajtmp != i){ 12084c1414c8SBarry Smith rs += PetscAbsScalar(*rtmp1++); 12094c1414c8SBarry Smith } else { 12104c1414c8SBarry Smith rs -= PetscRealPart(*rtmp1++); 12114c1414c8SBarry Smith } 12124c1414c8SBarry Smith ajtmp++; 12134c1414c8SBarry Smith } 12144c1414c8SBarry Smith if (rs>sctx.shift_top) sctx.shift_top = rs; 12154c1414c8SBarry Smith } 12164c1414c8SBarry Smith if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12; 12174c1414c8SBarry Smith sctx.shift_top *= 1.1; 12184c1414c8SBarry Smith sctx.nshift_max = 5; 12194c1414c8SBarry Smith sctx.shift_lo = 0.; 12204c1414c8SBarry Smith sctx.shift_hi = 1.; 12214c1414c8SBarry Smith } 12224c1414c8SBarry Smith sctx.shift_amount = 0; 12234c1414c8SBarry Smith sctx.nshift = 0; 12244c1414c8SBarry Smith 12254c1414c8SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 12264c1414c8SBarry Smith ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 12274c1414c8SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1228d9fead3dSBarry Smith ierr = PetscMalloc((3*n+1)*sizeof(PetscScalar),&rtmp11);CHKERRQ(ierr); 1229d9fead3dSBarry Smith ierr = PetscMemzero(rtmp11,(3*n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 12304c1414c8SBarry Smith ics = ic ; 1231d9fead3dSBarry Smith rtmp22 = rtmp11 + n; 1232d9fead3dSBarry Smith rtmp33 = rtmp22 + n; 12334c1414c8SBarry Smith 12344c1414c8SBarry Smith node_max = a->inode.node_count; 12354c1414c8SBarry Smith ns = a->inode.size; 12364c1414c8SBarry Smith if (!ns){ 12374c1414c8SBarry Smith SETERRQ(PETSC_ERR_PLIB,"Matrix without inode information"); 12384c1414c8SBarry Smith } 12394c1414c8SBarry Smith 12404c1414c8SBarry Smith /* If max inode size > 3, split it into two inodes.*/ 12414c1414c8SBarry Smith /* also map the inode sizes according to the ordering */ 12424c1414c8SBarry Smith ierr = PetscMalloc((n+1)* sizeof(PetscInt),&tmp_vec1);CHKERRQ(ierr); 12434c1414c8SBarry Smith for (i=0,j=0; i<node_max; ++i,++j){ 12444c1414c8SBarry Smith if (ns[i]>3) { 12454c1414c8SBarry Smith tmp_vec1[j] = ns[i]/2; /* Assuming ns[i] < =5 */ 12464c1414c8SBarry Smith ++j; 12474c1414c8SBarry Smith tmp_vec1[j] = ns[i] - tmp_vec1[j-1]; 12484c1414c8SBarry Smith } else { 12494c1414c8SBarry Smith tmp_vec1[j] = ns[i]; 12504c1414c8SBarry Smith } 12514c1414c8SBarry Smith } 12524c1414c8SBarry Smith /* Use the correct node_max */ 12534c1414c8SBarry Smith node_max = j; 12544c1414c8SBarry Smith 12554c1414c8SBarry Smith /* Now reorder the inode info based on mat re-ordering info */ 12564c1414c8SBarry Smith /* First create a row -> inode_size_array_index map */ 12574c1414c8SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt)+1,&nsmap);CHKERRQ(ierr); 12584c1414c8SBarry Smith ierr = PetscMalloc(node_max*sizeof(PetscInt)+1,&tmp_vec2);CHKERRQ(ierr); 12594c1414c8SBarry Smith for (i=0,row=0; i<node_max; i++) { 12604c1414c8SBarry Smith nodesz = tmp_vec1[i]; 12614c1414c8SBarry Smith for (j=0; j<nodesz; j++,row++) { 12624c1414c8SBarry Smith nsmap[row] = i; 12634c1414c8SBarry Smith } 12644c1414c8SBarry Smith } 12654c1414c8SBarry Smith /* Using nsmap, create a reordered ns structure */ 12664c1414c8SBarry Smith for (i=0,j=0; i< node_max; i++) { 12674c1414c8SBarry Smith nodesz = tmp_vec1[nsmap[r[j]]]; /* here the reordered row_no is in r[] */ 12684c1414c8SBarry Smith tmp_vec2[i] = nodesz; 12694c1414c8SBarry Smith j += nodesz; 12704c1414c8SBarry Smith } 12714c1414c8SBarry Smith ierr = PetscFree(nsmap);CHKERRQ(ierr); 12724c1414c8SBarry Smith ierr = PetscFree(tmp_vec1);CHKERRQ(ierr); 12734c1414c8SBarry Smith /* Now use the correct ns */ 12744c1414c8SBarry Smith ns = tmp_vec2; 12754c1414c8SBarry Smith 12764c1414c8SBarry Smith do { 12774c1414c8SBarry Smith sctx.lushift = PETSC_FALSE; 12784c1414c8SBarry Smith /* Now loop over each block-row, and do the factorization */ 12794c1414c8SBarry Smith for (i=0,row=0; i<node_max; i++) { 12804c1414c8SBarry Smith nodesz = ns[i]; 12814c1414c8SBarry Smith nz = bi[row+1] - bi[row]; 12824c1414c8SBarry Smith bjtmp = bj + bi[row]; 12834c1414c8SBarry Smith 12844c1414c8SBarry Smith switch (nodesz){ 12854c1414c8SBarry Smith case 1: 12864c1414c8SBarry Smith for (j=0; j<nz; j++){ 12874c1414c8SBarry Smith idx = bjtmp[j]; 1288d9fead3dSBarry Smith rtmp11[idx] = 0.0; 12894c1414c8SBarry Smith } 12904c1414c8SBarry Smith 12914c1414c8SBarry Smith /* load in initial (unfactored row) */ 12924c1414c8SBarry Smith idx = r[row]; 12934c1414c8SBarry Smith nz_tmp = ai[idx+1] - ai[idx]; 12944c1414c8SBarry Smith ajtmp = aj + ai[idx]; 12954c1414c8SBarry Smith v1 = aa + ai[idx]; 12964c1414c8SBarry Smith 12974c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 12984c1414c8SBarry Smith idx = ics[ajtmp[j]]; 1299d9fead3dSBarry Smith rtmp11[idx] = v1[j]; 13004c1414c8SBarry Smith } 1301d9fead3dSBarry Smith rtmp11[ics[r[row]]] += sctx.shift_amount; 13024c1414c8SBarry Smith 13034c1414c8SBarry Smith prow = *bjtmp++ ; 13044c1414c8SBarry Smith while (prow < row) { 1305d9fead3dSBarry Smith pc1 = rtmp11 + prow; 13064c1414c8SBarry Smith if (*pc1 != 0.0){ 13074c1414c8SBarry Smith pv = ba + bd[prow]; 13084c1414c8SBarry Smith pj = nbj + bd[prow]; 13094c1414c8SBarry Smith mul1 = *pc1 * *pv++; 13104c1414c8SBarry Smith *pc1 = mul1; 13114c1414c8SBarry Smith nz_tmp = bi[prow+1] - bd[prow] - 1; 13124c1414c8SBarry Smith ierr = PetscLogFlops(2*nz_tmp);CHKERRQ(ierr); 13134c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 13144c1414c8SBarry Smith tmp = pv[j]; 13154c1414c8SBarry Smith idx = pj[j]; 1316d9fead3dSBarry Smith rtmp11[idx] -= mul1 * tmp; 13174c1414c8SBarry Smith } 13184c1414c8SBarry Smith } 13194c1414c8SBarry Smith prow = *bjtmp++ ; 13204c1414c8SBarry Smith } 13214c1414c8SBarry Smith pj = bj + bi[row]; 13224c1414c8SBarry Smith pc1 = ba + bi[row]; 13234c1414c8SBarry Smith 1324d9fead3dSBarry Smith sctx.pv = rtmp11[row]; 1325d9fead3dSBarry Smith rtmp11[row] = 1.0/rtmp11[row]; /* invert diag */ 13264c1414c8SBarry Smith rs = 0.0; 13274c1414c8SBarry Smith for (j=0; j<nz; j++) { 13284c1414c8SBarry Smith idx = pj[j]; 1329d9fead3dSBarry Smith pc1[j] = rtmp11[idx]; /* rtmp11 -> ba */ 13304c1414c8SBarry Smith if (idx != row) rs += PetscAbsScalar(pc1[j]); 13314c1414c8SBarry Smith } 13324c1414c8SBarry Smith sctx.rs = rs; 1333c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,row,newshift);CHKERRQ(ierr); 13344c1414c8SBarry Smith if (newshift == 1) goto endofwhile; 13354c1414c8SBarry Smith break; 13364c1414c8SBarry Smith 13374c1414c8SBarry Smith case 2: 13384c1414c8SBarry Smith for (j=0; j<nz; j++) { 13394c1414c8SBarry Smith idx = bjtmp[j]; 1340d9fead3dSBarry Smith rtmp11[idx] = 0.0; 1341d9fead3dSBarry Smith rtmp22[idx] = 0.0; 13424c1414c8SBarry Smith } 13434c1414c8SBarry Smith 13444c1414c8SBarry Smith /* load in initial (unfactored row) */ 13454c1414c8SBarry Smith idx = r[row]; 13464c1414c8SBarry Smith nz_tmp = ai[idx+1] - ai[idx]; 13474c1414c8SBarry Smith ajtmp = aj + ai[idx]; 13484c1414c8SBarry Smith v1 = aa + ai[idx]; 13494c1414c8SBarry Smith v2 = aa + ai[idx+1]; 13504c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 13514c1414c8SBarry Smith idx = ics[ajtmp[j]]; 1352d9fead3dSBarry Smith rtmp11[idx] = v1[j]; 1353d9fead3dSBarry Smith rtmp22[idx] = v2[j]; 13544c1414c8SBarry Smith } 1355d9fead3dSBarry Smith rtmp11[ics[r[row]]] += sctx.shift_amount; 1356d9fead3dSBarry Smith rtmp22[ics[r[row+1]]] += sctx.shift_amount; 13574c1414c8SBarry Smith 13584c1414c8SBarry Smith prow = *bjtmp++ ; 13594c1414c8SBarry Smith while (prow < row) { 1360d9fead3dSBarry Smith pc1 = rtmp11 + prow; 1361d9fead3dSBarry Smith pc2 = rtmp22 + prow; 13624c1414c8SBarry Smith if (*pc1 != 0.0 || *pc2 != 0.0){ 13634c1414c8SBarry Smith pv = ba + bd[prow]; 13644c1414c8SBarry Smith pj = nbj + bd[prow]; 13654c1414c8SBarry Smith mul1 = *pc1 * *pv; 13664c1414c8SBarry Smith mul2 = *pc2 * *pv; 13674c1414c8SBarry Smith ++pv; 13684c1414c8SBarry Smith *pc1 = mul1; 13694c1414c8SBarry Smith *pc2 = mul2; 13704c1414c8SBarry Smith 13714c1414c8SBarry Smith nz_tmp = bi[prow+1] - bd[prow] - 1; 13724c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 13734c1414c8SBarry Smith tmp = pv[j]; 13744c1414c8SBarry Smith idx = pj[j]; 1375d9fead3dSBarry Smith rtmp11[idx] -= mul1 * tmp; 1376d9fead3dSBarry Smith rtmp22[idx] -= mul2 * tmp; 13774c1414c8SBarry Smith } 13784c1414c8SBarry Smith ierr = PetscLogFlops(4*nz_tmp);CHKERRQ(ierr); 13794c1414c8SBarry Smith } 13804c1414c8SBarry Smith prow = *bjtmp++ ; 13814c1414c8SBarry Smith } 13824c1414c8SBarry Smith 13834c1414c8SBarry Smith /* Now take care of diagonal 2x2 block. Note: prow = row here */ 1384d9fead3dSBarry Smith pc1 = rtmp11 + prow; 1385d9fead3dSBarry Smith pc2 = rtmp22 + prow; 13864c1414c8SBarry Smith 13874c1414c8SBarry Smith sctx.pv = *pc1; 13884c1414c8SBarry Smith pj = bj + bi[prow]; 13894c1414c8SBarry Smith rs = 0.0; 13904c1414c8SBarry Smith for (j=0; j<nz; j++){ 13914c1414c8SBarry Smith idx = pj[j]; 1392d9fead3dSBarry Smith if (idx != prow) rs += PetscAbsScalar(rtmp11[idx]); 13934c1414c8SBarry Smith } 13944c1414c8SBarry Smith sctx.rs = rs; 1395c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,row,newshift);CHKERRQ(ierr); 13964c1414c8SBarry Smith if (newshift == 1) goto endofwhile; 13974c1414c8SBarry Smith 13984c1414c8SBarry Smith if (*pc2 != 0.0){ 13994c1414c8SBarry Smith pj = nbj + bd[prow]; 14004c1414c8SBarry Smith mul2 = (*pc2)/(*pc1); /* since diag is not yet inverted.*/ 14014c1414c8SBarry Smith *pc2 = mul2; 14024c1414c8SBarry Smith nz_tmp = bi[prow+1] - bd[prow] - 1; 14034c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 14044c1414c8SBarry Smith idx = pj[j] ; 1405d9fead3dSBarry Smith tmp = rtmp11[idx]; 1406d9fead3dSBarry Smith rtmp22[idx] -= mul2 * tmp; 14074c1414c8SBarry Smith } 14084c1414c8SBarry Smith ierr = PetscLogFlops(2*nz_tmp);CHKERRQ(ierr); 14094c1414c8SBarry Smith } 14104c1414c8SBarry Smith 14114c1414c8SBarry Smith pj = bj + bi[row]; 14124c1414c8SBarry Smith pc1 = ba + bi[row]; 14134c1414c8SBarry Smith pc2 = ba + bi[row+1]; 14144c1414c8SBarry Smith 1415d9fead3dSBarry Smith sctx.pv = rtmp22[row+1]; 14164c1414c8SBarry Smith rs = 0.0; 1417d9fead3dSBarry Smith rtmp11[row] = 1.0/rtmp11[row]; 1418d9fead3dSBarry Smith rtmp22[row+1] = 1.0/rtmp22[row+1]; 14194c1414c8SBarry Smith /* copy row entries from dense representation to sparse */ 14204c1414c8SBarry Smith for (j=0; j<nz; j++) { 14214c1414c8SBarry Smith idx = pj[j]; 1422d9fead3dSBarry Smith pc1[j] = rtmp11[idx]; 1423d9fead3dSBarry Smith pc2[j] = rtmp22[idx]; 14244c1414c8SBarry Smith if (idx != row+1) rs += PetscAbsScalar(pc2[j]); 14254c1414c8SBarry Smith } 14264c1414c8SBarry Smith sctx.rs = rs; 1427c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,row+1,newshift);CHKERRQ(ierr); 14284c1414c8SBarry Smith if (newshift == 1) goto endofwhile; 14294c1414c8SBarry Smith break; 14304c1414c8SBarry Smith 14314c1414c8SBarry Smith case 3: 14324c1414c8SBarry Smith for (j=0; j<nz; j++) { 14334c1414c8SBarry Smith idx = bjtmp[j]; 1434d9fead3dSBarry Smith rtmp11[idx] = 0.0; 1435d9fead3dSBarry Smith rtmp22[idx] = 0.0; 1436d9fead3dSBarry Smith rtmp33[idx] = 0.0; 14374c1414c8SBarry Smith } 14384c1414c8SBarry Smith /* copy the nonzeros for the 3 rows from sparse representation to dense in rtmp*[] */ 14394c1414c8SBarry Smith idx = r[row]; 14404c1414c8SBarry Smith nz_tmp = ai[idx+1] - ai[idx]; 14414c1414c8SBarry Smith ajtmp = aj + ai[idx]; 14424c1414c8SBarry Smith v1 = aa + ai[idx]; 14434c1414c8SBarry Smith v2 = aa + ai[idx+1]; 14444c1414c8SBarry Smith v3 = aa + ai[idx+2]; 14454c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 14464c1414c8SBarry Smith idx = ics[ajtmp[j]]; 1447d9fead3dSBarry Smith rtmp11[idx] = v1[j]; 1448d9fead3dSBarry Smith rtmp22[idx] = v2[j]; 1449d9fead3dSBarry Smith rtmp33[idx] = v3[j]; 14504c1414c8SBarry Smith } 1451d9fead3dSBarry Smith rtmp11[ics[r[row]]] += sctx.shift_amount; 1452d9fead3dSBarry Smith rtmp22[ics[r[row+1]]] += sctx.shift_amount; 1453d9fead3dSBarry Smith rtmp33[ics[r[row+2]]] += sctx.shift_amount; 14544c1414c8SBarry Smith 14554c1414c8SBarry Smith /* loop over all pivot row blocks above this row block */ 14564c1414c8SBarry Smith prow = *bjtmp++ ; 14574c1414c8SBarry Smith while (prow < row) { 1458d9fead3dSBarry Smith pc1 = rtmp11 + prow; 1459d9fead3dSBarry Smith pc2 = rtmp22 + prow; 1460d9fead3dSBarry Smith pc3 = rtmp33 + prow; 14614c1414c8SBarry Smith if (*pc1 != 0.0 || *pc2 != 0.0 || *pc3 !=0.0){ 14624c1414c8SBarry Smith pv = ba + bd[prow]; 14634c1414c8SBarry Smith pj = nbj + bd[prow]; 14644c1414c8SBarry Smith mul1 = *pc1 * *pv; 14654c1414c8SBarry Smith mul2 = *pc2 * *pv; 14664c1414c8SBarry Smith mul3 = *pc3 * *pv; 14674c1414c8SBarry Smith ++pv; 14684c1414c8SBarry Smith *pc1 = mul1; 14694c1414c8SBarry Smith *pc2 = mul2; 14704c1414c8SBarry Smith *pc3 = mul3; 14714c1414c8SBarry Smith 14724c1414c8SBarry Smith nz_tmp = bi[prow+1] - bd[prow] - 1; 14734c1414c8SBarry Smith /* update this row based on pivot row */ 14744c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 14754c1414c8SBarry Smith tmp = pv[j]; 14764c1414c8SBarry Smith idx = pj[j]; 1477d9fead3dSBarry Smith rtmp11[idx] -= mul1 * tmp; 1478d9fead3dSBarry Smith rtmp22[idx] -= mul2 * tmp; 1479d9fead3dSBarry Smith rtmp33[idx] -= mul3 * tmp; 14804c1414c8SBarry Smith } 14814c1414c8SBarry Smith ierr = PetscLogFlops(6*nz_tmp);CHKERRQ(ierr); 14824c1414c8SBarry Smith } 14834c1414c8SBarry Smith prow = *bjtmp++ ; 14844c1414c8SBarry Smith } 14854c1414c8SBarry Smith 14864c1414c8SBarry Smith /* Now take care of diagonal 3x3 block in this set of rows */ 14874c1414c8SBarry Smith /* note: prow = row here */ 1488d9fead3dSBarry Smith pc1 = rtmp11 + prow; 1489d9fead3dSBarry Smith pc2 = rtmp22 + prow; 1490d9fead3dSBarry Smith pc3 = rtmp33 + prow; 14914c1414c8SBarry Smith 14924c1414c8SBarry Smith sctx.pv = *pc1; 14934c1414c8SBarry Smith pj = bj + bi[prow]; 14944c1414c8SBarry Smith rs = 0.0; 14954c1414c8SBarry Smith for (j=0; j<nz; j++){ 14964c1414c8SBarry Smith idx = pj[j]; 1497d9fead3dSBarry Smith if (idx != row) rs += PetscAbsScalar(rtmp11[idx]); 14984c1414c8SBarry Smith } 14994c1414c8SBarry Smith sctx.rs = rs; 1500c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,row,newshift);CHKERRQ(ierr); 15014c1414c8SBarry Smith if (newshift == 1) goto endofwhile; 15024c1414c8SBarry Smith 15034c1414c8SBarry Smith if (*pc2 != 0.0 || *pc3 != 0.0){ 15044c1414c8SBarry Smith mul2 = (*pc2)/(*pc1); 15054c1414c8SBarry Smith mul3 = (*pc3)/(*pc1); 15064c1414c8SBarry Smith *pc2 = mul2; 15074c1414c8SBarry Smith *pc3 = mul3; 15084c1414c8SBarry Smith nz_tmp = bi[prow+1] - bd[prow] - 1; 15094c1414c8SBarry Smith pj = nbj + bd[prow]; 15104c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 15114c1414c8SBarry Smith idx = pj[j] ; 1512d9fead3dSBarry Smith tmp = rtmp11[idx]; 1513d9fead3dSBarry Smith rtmp22[idx] -= mul2 * tmp; 1514d9fead3dSBarry Smith rtmp33[idx] -= mul3 * tmp; 15154c1414c8SBarry Smith } 15164c1414c8SBarry Smith ierr = PetscLogFlops(4*nz_tmp);CHKERRQ(ierr); 15174c1414c8SBarry Smith } 15184c1414c8SBarry Smith ++prow; 15194c1414c8SBarry Smith 1520d9fead3dSBarry Smith pc2 = rtmp22 + prow; 1521d9fead3dSBarry Smith pc3 = rtmp33 + prow; 15224c1414c8SBarry Smith sctx.pv = *pc2; 15234c1414c8SBarry Smith pj = bj + bi[prow]; 15244c1414c8SBarry Smith rs = 0.0; 15254c1414c8SBarry Smith for (j=0; j<nz; j++){ 15264c1414c8SBarry Smith idx = pj[j]; 1527d9fead3dSBarry Smith if (idx != prow) rs += PetscAbsScalar(rtmp22[idx]); 15284c1414c8SBarry Smith } 15294c1414c8SBarry Smith sctx.rs = rs; 1530c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,row+1,newshift);CHKERRQ(ierr); 15314c1414c8SBarry Smith if (newshift == 1) goto endofwhile; 15324c1414c8SBarry Smith 15334c1414c8SBarry Smith if (*pc3 != 0.0){ 15344c1414c8SBarry Smith mul3 = (*pc3)/(*pc2); 15354c1414c8SBarry Smith *pc3 = mul3; 15364c1414c8SBarry Smith pj = nbj + bd[prow]; 15374c1414c8SBarry Smith nz_tmp = bi[prow+1] - bd[prow] - 1; 15384c1414c8SBarry Smith for (j=0; j<nz_tmp; j++) { 15394c1414c8SBarry Smith idx = pj[j] ; 1540d9fead3dSBarry Smith tmp = rtmp22[idx]; 1541d9fead3dSBarry Smith rtmp33[idx] -= mul3 * tmp; 15424c1414c8SBarry Smith } 15434c1414c8SBarry Smith ierr = PetscLogFlops(4*nz_tmp);CHKERRQ(ierr); 15444c1414c8SBarry Smith } 15454c1414c8SBarry Smith 15464c1414c8SBarry Smith pj = bj + bi[row]; 15474c1414c8SBarry Smith pc1 = ba + bi[row]; 15484c1414c8SBarry Smith pc2 = ba + bi[row+1]; 15494c1414c8SBarry Smith pc3 = ba + bi[row+2]; 15504c1414c8SBarry Smith 1551d9fead3dSBarry Smith sctx.pv = rtmp33[row+2]; 15524c1414c8SBarry Smith rs = 0.0; 1553d9fead3dSBarry Smith rtmp11[row] = 1.0/rtmp11[row]; 1554d9fead3dSBarry Smith rtmp22[row+1] = 1.0/rtmp22[row+1]; 1555d9fead3dSBarry Smith rtmp33[row+2] = 1.0/rtmp33[row+2]; 15564c1414c8SBarry Smith /* copy row entries from dense representation to sparse */ 15574c1414c8SBarry Smith for (j=0; j<nz; j++) { 15584c1414c8SBarry Smith idx = pj[j]; 1559d9fead3dSBarry Smith pc1[j] = rtmp11[idx]; 1560d9fead3dSBarry Smith pc2[j] = rtmp22[idx]; 1561d9fead3dSBarry Smith pc3[j] = rtmp33[idx]; 15624c1414c8SBarry Smith if (idx != row+2) rs += PetscAbsScalar(pc3[j]); 15634c1414c8SBarry Smith } 15644c1414c8SBarry Smith 15654c1414c8SBarry Smith sctx.rs = rs; 1566c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,row+2,newshift);CHKERRQ(ierr); 15674c1414c8SBarry Smith if (newshift == 1) goto endofwhile; 15684c1414c8SBarry Smith break; 15694c1414c8SBarry Smith 15704c1414c8SBarry Smith default: 15714c1414c8SBarry Smith SETERRQ(PETSC_ERR_COR,"Node size not yet supported \n"); 15724c1414c8SBarry Smith } 15734c1414c8SBarry Smith row += nodesz; /* Update the row */ 15744c1414c8SBarry Smith } 15754c1414c8SBarry Smith endofwhile:; 15764c1414c8SBarry Smith } while (sctx.lushift); 1577d9fead3dSBarry Smith ierr = PetscFree(rtmp11);CHKERRQ(ierr); 15784c1414c8SBarry Smith ierr = PetscFree(tmp_vec2);CHKERRQ(ierr); 15794c1414c8SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 15804c1414c8SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 15814c1414c8SBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 1582*719d5645SBarry Smith (B)->ops->solve = MatSolve_Inode; 1583db4efbfdSBarry Smith /* do not set solve add, since MatSolve_Inode + Add is faster */ 1584db4efbfdSBarry Smith C->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 1585db4efbfdSBarry Smith C->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ; 15864c1414c8SBarry Smith C->assembled = PETSC_TRUE; 15875c9eb25fSBarry Smith C->preallocated = PETSC_TRUE; 15884c1414c8SBarry Smith if (sctx.nshift) { 15894c1414c8SBarry Smith if (info->shiftnz) { 15901e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 15914c1414c8SBarry Smith } else if (info->shiftpd) { 15921e2582c4SBarry Smith ierr = PetscInfo4(A,"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); 15934c1414c8SBarry Smith } 15944c1414c8SBarry Smith } 1595d0f46423SBarry Smith ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 15964c1414c8SBarry Smith PetscFunctionReturn(0); 15974c1414c8SBarry Smith } 15984c1414c8SBarry Smith 15994c1414c8SBarry Smith /* 16004c1414c8SBarry Smith Makes a longer coloring[] array and calls the usual code with that 16014c1414c8SBarry Smith */ 16024c1414c8SBarry Smith #undef __FUNCT__ 16034c1414c8SBarry Smith #define __FUNCT__ "MatColoringPatch_Inode" 16044c1414c8SBarry Smith PetscErrorCode MatColoringPatch_Inode(Mat mat,PetscInt ncolors,PetscInt nin,ISColoringValue coloring[],ISColoring *iscoloring) 16054c1414c8SBarry Smith { 16064c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)mat->data; 16074c1414c8SBarry Smith PetscErrorCode ierr; 1608d0f46423SBarry Smith PetscInt n = mat->cmap->n,m = a->inode.node_count,j,*ns = a->inode.size,row; 16094c1414c8SBarry Smith PetscInt *colorused,i; 16104c1414c8SBarry Smith ISColoringValue *newcolor; 16114c1414c8SBarry Smith 16124c1414c8SBarry Smith PetscFunctionBegin; 16134c1414c8SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&newcolor);CHKERRQ(ierr); 16144c1414c8SBarry Smith /* loop over inodes, marking a color for each column*/ 16154c1414c8SBarry Smith row = 0; 16164c1414c8SBarry Smith for (i=0; i<m; i++){ 16174c1414c8SBarry Smith for (j=0; j<ns[i]; j++) { 16184c1414c8SBarry Smith newcolor[row++] = coloring[i] + j*ncolors; 16194c1414c8SBarry Smith } 16204c1414c8SBarry Smith } 16214c1414c8SBarry Smith 16224c1414c8SBarry Smith /* eliminate unneeded colors */ 16234c1414c8SBarry Smith ierr = PetscMalloc(5*ncolors*sizeof(PetscInt),&colorused);CHKERRQ(ierr); 16244c1414c8SBarry Smith ierr = PetscMemzero(colorused,5*ncolors*sizeof(PetscInt));CHKERRQ(ierr); 16254c1414c8SBarry Smith for (i=0; i<n; i++) { 16264c1414c8SBarry Smith colorused[newcolor[i]] = 1; 16274c1414c8SBarry Smith } 16284c1414c8SBarry Smith 16294c1414c8SBarry Smith for (i=1; i<5*ncolors; i++) { 16304c1414c8SBarry Smith colorused[i] += colorused[i-1]; 16314c1414c8SBarry Smith } 16324c1414c8SBarry Smith ncolors = colorused[5*ncolors-1]; 16334c1414c8SBarry Smith for (i=0; i<n; i++) { 163457815681SSatish Balay newcolor[i] = colorused[newcolor[i]]-1; 16354c1414c8SBarry Smith } 16364c1414c8SBarry Smith ierr = PetscFree(colorused);CHKERRQ(ierr); 16377adad957SLisandro Dalcin ierr = ISColoringCreate(((PetscObject)mat)->comm,ncolors,n,newcolor,iscoloring);CHKERRQ(ierr); 16384c1414c8SBarry Smith ierr = PetscFree(coloring);CHKERRQ(ierr); 16394c1414c8SBarry Smith PetscFunctionReturn(0); 16404c1414c8SBarry Smith } 16414c1414c8SBarry Smith 16422af78befSBarry Smith #include "src/inline/ilu.h" 16432af78befSBarry Smith 16442af78befSBarry Smith #undef __FUNCT__ 16452af78befSBarry Smith #define __FUNCT__ "MatRelax_Inode" 16462af78befSBarry Smith PetscErrorCode MatRelax_Inode(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 16472af78befSBarry Smith { 16482af78befSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1649a77337e4SBarry Smith PetscScalar *x,*xs,sum1,sum2,sum3,sum4,sum5,tmp0,tmp1,tmp2,tmp3; 1650a77337e4SBarry Smith MatScalar *ibdiag,*bdiag; 1651d9fead3dSBarry Smith PetscScalar *b,*xb,tmp4,tmp5,x1,x2,x3,x4,x5; 1652dd6ea824SBarry Smith const MatScalar *v = a->a,*v1,*v2,*v3,*v4,*v5; 165362bba022SBarry Smith PetscReal zeropivot = 1.0e-15, shift = 0.0; 16542af78befSBarry Smith PetscErrorCode ierr; 1655f0d39aaaSBarry Smith PetscInt n,m = a->inode.node_count,*sizes = a->inode.size,cnt = 0,i,j,row,i1,i2; 1656f0d39aaaSBarry Smith PetscInt *idx,*diag = a->diag,*ii = a->i,sz,k; 16572af78befSBarry Smith 16582af78befSBarry Smith PetscFunctionBegin; 165971f1c65dSBarry Smith if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"No support for omega != 1.0; use -mat_no_inode"); 166071f1c65dSBarry Smith if (fshift != 0.0) SETERRQ(PETSC_ERR_SUP,"No support for fshift != 0.0; use -mat_no_inode"); 166171f1c65dSBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support for Eisenstat trick; use -mat_no_inode"); 16622af78befSBarry Smith 166371f1c65dSBarry Smith if (!a->inode.ibdiagvalid) { 16642af78befSBarry Smith if (!a->inode.ibdiag) { 16652af78befSBarry Smith /* calculate space needed for diagonal blocks */ 16662af78befSBarry Smith for (i=0; i<m; i++) { 16672af78befSBarry Smith cnt += sizes[i]*sizes[i]; 16682af78befSBarry Smith } 1669f0d39aaaSBarry Smith a->inode.bdiagsize = cnt; 1670a77337e4SBarry Smith ierr = PetscMalloc2(cnt,MatScalar,&a->inode.ibdiag,cnt,MatScalar,&a->inode.bdiag);CHKERRQ(ierr); 167171f1c65dSBarry Smith } 167271f1c65dSBarry Smith 167371f1c65dSBarry Smith /* copy over the diagonal blocks and invert them */ 16742af78befSBarry Smith ibdiag = a->inode.ibdiag; 16752af78befSBarry Smith bdiag = a->inode.bdiag; 16762af78befSBarry Smith cnt = 0; 16772af78befSBarry Smith for (i=0, row = 0; i<m; i++) { 16782af78befSBarry Smith for (j=0; j<sizes[i]; j++) { 1679f0d39aaaSBarry Smith for (k=0; k<sizes[i]; k++) { 1680f0d39aaaSBarry Smith bdiag[cnt+k*sizes[i]+j] = v[diag[row+j] - j + k]; 1681f0d39aaaSBarry Smith } 16822af78befSBarry Smith } 1683a77337e4SBarry Smith ierr = PetscMemcpy(ibdiag+cnt,bdiag+cnt,sizes[i]*sizes[i]*sizeof(MatScalar));CHKERRQ(ierr); 16842af78befSBarry Smith 16852af78befSBarry Smith switch(sizes[i]) { 16862af78befSBarry Smith case 1: 16872af78befSBarry Smith /* Create matrix data structure */ 168864c62002SMatthew Knepley if (PetscAbsScalar(ibdiag[cnt]) < zeropivot) SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot on row %D",row); 168964c62002SMatthew Knepley ibdiag[cnt] = 1.0/ibdiag[cnt]; 16902af78befSBarry Smith break; 16912af78befSBarry Smith case 2: 169262bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(ibdiag+cnt,shift);CHKERRQ(ierr); 16932af78befSBarry Smith break; 16942af78befSBarry Smith case 3: 169562bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_3(ibdiag+cnt,shift);CHKERRQ(ierr); 16962af78befSBarry Smith break; 16972af78befSBarry Smith case 4: 169862bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_4(ibdiag+cnt,shift);CHKERRQ(ierr); 16992af78befSBarry Smith break; 17002af78befSBarry Smith case 5: 170162bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_5(ibdiag+cnt,shift);CHKERRQ(ierr); 17022af78befSBarry Smith break; 17032af78befSBarry Smith default: 17042af78befSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Inode size %D not supported",sizes[i]); 17052af78befSBarry Smith } 17062af78befSBarry Smith cnt += sizes[i]*sizes[i]; 17072af78befSBarry Smith row += sizes[i]; 17082af78befSBarry Smith } 170971f1c65dSBarry Smith a->inode.ibdiagvalid = PETSC_TRUE; 17102af78befSBarry Smith } 17112af78befSBarry Smith ibdiag = a->inode.ibdiag; 17122af78befSBarry Smith bdiag = a->inode.bdiag; 17132af78befSBarry Smith 17142af78befSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 17152af78befSBarry Smith if (xx != bb) { 17162af78befSBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 17172af78befSBarry Smith } else { 17182af78befSBarry Smith b = x; 17192af78befSBarry Smith } 17202af78befSBarry Smith 17212af78befSBarry Smith /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */ 17222af78befSBarry Smith xs = x; 17232af78befSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 17242af78befSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 17252af78befSBarry Smith 17268862d2efSBarry Smith for (i=0, row=0; i<m; i++) { 17278862d2efSBarry Smith sz = diag[row] - ii[row]; 17288862d2efSBarry Smith v1 = a->a + ii[row]; 17298862d2efSBarry Smith idx = a->j + ii[row]; 17308862d2efSBarry Smith 17318862d2efSBarry Smith /* see comments for MatMult_Inode() for how this is coded */ 17328862d2efSBarry Smith switch (sizes[i]){ 17338862d2efSBarry Smith case 1: 17348862d2efSBarry Smith 17358862d2efSBarry Smith sum1 = b[row]; 17368862d2efSBarry Smith for(n = 0; n<sz-1; n+=2) { 17378862d2efSBarry Smith i1 = idx[0]; 17388862d2efSBarry Smith i2 = idx[1]; 17398862d2efSBarry Smith idx += 2; 17408862d2efSBarry Smith tmp0 = x[i1]; 17418862d2efSBarry Smith tmp1 = x[i2]; 17428862d2efSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 17438862d2efSBarry Smith } 17448862d2efSBarry Smith 17458862d2efSBarry Smith if (n == sz-1){ 1746f0d39aaaSBarry Smith tmp0 = x[*idx]; 1747f0d39aaaSBarry Smith sum1 -= *v1 * tmp0; 17488862d2efSBarry Smith } 17498862d2efSBarry Smith x[row++] = sum1*(*ibdiag++); 17508862d2efSBarry Smith break; 1751f0d39aaaSBarry Smith case 2: 1752f0d39aaaSBarry Smith v2 = a->a + ii[row+1]; 1753f0d39aaaSBarry Smith sum1 = b[row]; 1754f0d39aaaSBarry Smith sum2 = b[row+1]; 1755f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 1756f0d39aaaSBarry Smith i1 = idx[0]; 1757f0d39aaaSBarry Smith i2 = idx[1]; 1758f0d39aaaSBarry Smith idx += 2; 1759f0d39aaaSBarry Smith tmp0 = x[i1]; 1760f0d39aaaSBarry Smith tmp1 = x[i2]; 1761f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 1762f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 1763f0d39aaaSBarry Smith } 1764f0d39aaaSBarry Smith 1765f0d39aaaSBarry Smith if (n == sz-1){ 1766f0d39aaaSBarry Smith tmp0 = x[*idx]; 1767f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0; 1768f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0; 1769f0d39aaaSBarry Smith } 1770f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[0] + sum2*ibdiag[2]; 1771f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[1] + sum2*ibdiag[3]; 1772f0d39aaaSBarry Smith ibdiag += 4; 1773f0d39aaaSBarry Smith break; 1774f0d39aaaSBarry Smith case 3: 1775f0d39aaaSBarry Smith v2 = a->a + ii[row+1]; 1776f0d39aaaSBarry Smith v3 = a->a + ii[row+2]; 1777f0d39aaaSBarry Smith sum1 = b[row]; 1778f0d39aaaSBarry Smith sum2 = b[row+1]; 1779f0d39aaaSBarry Smith sum3 = b[row+2]; 1780f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 1781f0d39aaaSBarry Smith i1 = idx[0]; 1782f0d39aaaSBarry Smith i2 = idx[1]; 1783f0d39aaaSBarry Smith idx += 2; 1784f0d39aaaSBarry Smith tmp0 = x[i1]; 1785f0d39aaaSBarry Smith tmp1 = x[i2]; 1786f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 1787f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 1788f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 1789f0d39aaaSBarry Smith } 1790f0d39aaaSBarry Smith 1791f0d39aaaSBarry Smith if (n == sz-1){ 1792f0d39aaaSBarry Smith tmp0 = x[*idx]; 1793f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0; 1794f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0; 1795f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0; 1796f0d39aaaSBarry Smith } 1797f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[0] + sum2*ibdiag[3] + sum3*ibdiag[6]; 1798f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[1] + sum2*ibdiag[4] + sum3*ibdiag[7]; 1799f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[2] + sum2*ibdiag[5] + sum3*ibdiag[8]; 1800f0d39aaaSBarry Smith ibdiag += 9; 1801f0d39aaaSBarry Smith break; 1802f0d39aaaSBarry Smith case 4: 1803f0d39aaaSBarry Smith v2 = a->a + ii[row+1]; 1804f0d39aaaSBarry Smith v3 = a->a + ii[row+2]; 1805f0d39aaaSBarry Smith v4 = a->a + ii[row+3]; 1806f0d39aaaSBarry Smith sum1 = b[row]; 1807f0d39aaaSBarry Smith sum2 = b[row+1]; 1808f0d39aaaSBarry Smith sum3 = b[row+2]; 1809f0d39aaaSBarry Smith sum4 = b[row+3]; 1810f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 1811f0d39aaaSBarry Smith i1 = idx[0]; 1812f0d39aaaSBarry Smith i2 = idx[1]; 1813f0d39aaaSBarry Smith idx += 2; 1814f0d39aaaSBarry Smith tmp0 = x[i1]; 1815f0d39aaaSBarry Smith tmp1 = x[i2]; 1816f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 1817f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 1818f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 1819f0d39aaaSBarry Smith sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2; 1820f0d39aaaSBarry Smith } 1821f0d39aaaSBarry Smith 1822f0d39aaaSBarry Smith if (n == sz-1){ 1823f0d39aaaSBarry Smith tmp0 = x[*idx]; 1824f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0; 1825f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0; 1826f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0; 1827f0d39aaaSBarry Smith sum4 -= v4[0] * tmp0; 1828f0d39aaaSBarry Smith } 1829f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[0] + sum2*ibdiag[4] + sum3*ibdiag[8] + sum4*ibdiag[12]; 1830f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[1] + sum2*ibdiag[5] + sum3*ibdiag[9] + sum4*ibdiag[13]; 1831f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[2] + sum2*ibdiag[6] + sum3*ibdiag[10] + sum4*ibdiag[14]; 1832f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[3] + sum2*ibdiag[7] + sum3*ibdiag[11] + sum4*ibdiag[15]; 1833f0d39aaaSBarry Smith ibdiag += 16; 1834f0d39aaaSBarry Smith break; 1835f0d39aaaSBarry Smith case 5: 1836f0d39aaaSBarry Smith v2 = a->a + ii[row+1]; 1837f0d39aaaSBarry Smith v3 = a->a + ii[row+2]; 1838f0d39aaaSBarry Smith v4 = a->a + ii[row+3]; 1839f0d39aaaSBarry Smith v5 = a->a + ii[row+4]; 1840f0d39aaaSBarry Smith sum1 = b[row]; 1841f0d39aaaSBarry Smith sum2 = b[row+1]; 1842f0d39aaaSBarry Smith sum3 = b[row+2]; 1843f0d39aaaSBarry Smith sum4 = b[row+3]; 1844f0d39aaaSBarry Smith sum5 = b[row+4]; 1845f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 1846f0d39aaaSBarry Smith i1 = idx[0]; 1847f0d39aaaSBarry Smith i2 = idx[1]; 1848f0d39aaaSBarry Smith idx += 2; 1849f0d39aaaSBarry Smith tmp0 = x[i1]; 1850f0d39aaaSBarry Smith tmp1 = x[i2]; 1851f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 1852f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 1853f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 1854f0d39aaaSBarry Smith sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2; 1855f0d39aaaSBarry Smith sum5 -= v5[0] * tmp0 + v5[1] * tmp1; v5 += 2; 1856f0d39aaaSBarry Smith } 1857f0d39aaaSBarry Smith 1858f0d39aaaSBarry Smith if (n == sz-1){ 1859f0d39aaaSBarry Smith tmp0 = x[*idx]; 1860f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0; 1861f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0; 1862f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0; 1863f0d39aaaSBarry Smith sum4 -= v4[0] * tmp0; 1864f0d39aaaSBarry Smith sum5 -= v5[0] * tmp0; 1865f0d39aaaSBarry Smith } 1866f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[0] + sum2*ibdiag[5] + sum3*ibdiag[10] + sum4*ibdiag[15] + sum5*ibdiag[20]; 1867f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[1] + sum2*ibdiag[6] + sum3*ibdiag[11] + sum4*ibdiag[16] + sum5*ibdiag[21]; 1868f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[2] + sum2*ibdiag[7] + sum3*ibdiag[12] + sum4*ibdiag[17] + sum5*ibdiag[22]; 1869f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[3] + sum2*ibdiag[8] + sum3*ibdiag[13] + sum4*ibdiag[18] + sum5*ibdiag[23]; 1870f0d39aaaSBarry Smith x[row++] = sum1*ibdiag[4] + sum2*ibdiag[9] + sum3*ibdiag[14] + sum4*ibdiag[19] + sum5*ibdiag[24]; 1871f0d39aaaSBarry Smith ibdiag += 25; 1872f0d39aaaSBarry Smith break; 18738862d2efSBarry Smith default: 18748862d2efSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Inode size %D not supported",sizes[i]); 18758862d2efSBarry Smith } 18762af78befSBarry Smith } 18772af78befSBarry Smith 18782af78befSBarry Smith xb = x; 18792af78befSBarry Smith ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 18802af78befSBarry Smith } else xb = b; 18812af78befSBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 18822af78befSBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 18838862d2efSBarry Smith cnt = 0; 18848862d2efSBarry Smith for (i=0, row=0; i<m; i++) { 18852af78befSBarry Smith 18868862d2efSBarry Smith switch (sizes[i]){ 18878862d2efSBarry Smith case 1: 18888862d2efSBarry Smith x[row++] *= bdiag[cnt++]; 18898862d2efSBarry Smith break; 1890f0d39aaaSBarry Smith case 2: 1891f0d39aaaSBarry Smith x1 = x[row]; x2 = x[row+1]; 1892f0d39aaaSBarry Smith tmp1 = x1*bdiag[cnt] + x2*bdiag[cnt+2]; 1893f0d39aaaSBarry Smith tmp2 = x1*bdiag[cnt+1] + x2*bdiag[cnt+3]; 1894f0d39aaaSBarry Smith x[row++] = tmp1; 1895f0d39aaaSBarry Smith x[row++] = tmp2; 1896f0d39aaaSBarry Smith cnt += 4; 1897f0d39aaaSBarry Smith break; 1898f0d39aaaSBarry Smith case 3: 1899f0d39aaaSBarry Smith x1 = x[row]; x2 = x[row+1]; x3 = x[row+2]; 1900f0d39aaaSBarry Smith tmp1 = x1*bdiag[cnt] + x2*bdiag[cnt+3] + x3*bdiag[cnt+6]; 1901f0d39aaaSBarry Smith tmp2 = x1*bdiag[cnt+1] + x2*bdiag[cnt+4] + x3*bdiag[cnt+7]; 1902f0d39aaaSBarry Smith tmp3 = x1*bdiag[cnt+2] + x2*bdiag[cnt+5] + x3*bdiag[cnt+8]; 1903f0d39aaaSBarry Smith x[row++] = tmp1; 1904f0d39aaaSBarry Smith x[row++] = tmp2; 1905f0d39aaaSBarry Smith x[row++] = tmp3; 1906f0d39aaaSBarry Smith cnt += 9; 1907f0d39aaaSBarry Smith break; 1908f0d39aaaSBarry Smith case 4: 1909f0d39aaaSBarry Smith x1 = x[row]; x2 = x[row+1]; x3 = x[row+2]; x4 = x[row+3]; 1910f0d39aaaSBarry Smith tmp1 = x1*bdiag[cnt] + x2*bdiag[cnt+4] + x3*bdiag[cnt+8] + x4*bdiag[cnt+12]; 1911f0d39aaaSBarry Smith tmp2 = x1*bdiag[cnt+1] + x2*bdiag[cnt+5] + x3*bdiag[cnt+9] + x4*bdiag[cnt+13]; 1912f0d39aaaSBarry Smith tmp3 = x1*bdiag[cnt+2] + x2*bdiag[cnt+6] + x3*bdiag[cnt+10] + x4*bdiag[cnt+14]; 1913f0d39aaaSBarry Smith tmp4 = x1*bdiag[cnt+3] + x2*bdiag[cnt+7] + x3*bdiag[cnt+11] + x4*bdiag[cnt+15]; 1914f0d39aaaSBarry Smith x[row++] = tmp1; 1915f0d39aaaSBarry Smith x[row++] = tmp2; 1916f0d39aaaSBarry Smith x[row++] = tmp3; 1917f0d39aaaSBarry Smith x[row++] = tmp4; 1918f0d39aaaSBarry Smith cnt += 16; 1919f0d39aaaSBarry Smith break; 1920f0d39aaaSBarry Smith case 5: 1921f0d39aaaSBarry Smith x1 = x[row]; x2 = x[row+1]; x3 = x[row+2]; x4 = x[row+3]; x5 = x[row+4]; 1922f0d39aaaSBarry Smith tmp1 = x1*bdiag[cnt] + x2*bdiag[cnt+5] + x3*bdiag[cnt+10] + x4*bdiag[cnt+15] + x5*bdiag[cnt+20]; 1923f0d39aaaSBarry Smith tmp2 = x1*bdiag[cnt+1] + x2*bdiag[cnt+6] + x3*bdiag[cnt+11] + x4*bdiag[cnt+16] + x5*bdiag[cnt+21]; 1924f0d39aaaSBarry Smith tmp3 = x1*bdiag[cnt+2] + x2*bdiag[cnt+7] + x3*bdiag[cnt+12] + x4*bdiag[cnt+17] + x5*bdiag[cnt+22]; 1925f0d39aaaSBarry Smith tmp4 = x1*bdiag[cnt+3] + x2*bdiag[cnt+8] + x3*bdiag[cnt+13] + x4*bdiag[cnt+18] + x5*bdiag[cnt+23]; 1926f0d39aaaSBarry Smith tmp5 = x1*bdiag[cnt+4] + x2*bdiag[cnt+9] + x3*bdiag[cnt+14] + x4*bdiag[cnt+19] + x5*bdiag[cnt+24]; 1927f0d39aaaSBarry Smith x[row++] = tmp1; 1928f0d39aaaSBarry Smith x[row++] = tmp2; 1929f0d39aaaSBarry Smith x[row++] = tmp3; 1930f0d39aaaSBarry Smith x[row++] = tmp4; 1931f0d39aaaSBarry Smith x[row++] = tmp5; 1932f0d39aaaSBarry Smith cnt += 25; 1933f0d39aaaSBarry Smith break; 19348862d2efSBarry Smith default: 19358862d2efSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Inode size %D not supported",sizes[i]); 19368862d2efSBarry Smith } 19372af78befSBarry Smith } 19382af78befSBarry Smith ierr = PetscLogFlops(m);CHKERRQ(ierr); 19392af78befSBarry Smith } 19402af78befSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 19412af78befSBarry Smith 1942f0d39aaaSBarry Smith ibdiag = a->inode.ibdiag+a->inode.bdiagsize; 1943d0f46423SBarry Smith for (i=m-1, row=A->rmap->n-1; i>=0; i--) { 1944f0d39aaaSBarry Smith ibdiag -= sizes[i]*sizes[i]; 19458862d2efSBarry Smith sz = ii[row+1] - diag[row] - 1; 19468862d2efSBarry Smith v1 = a->a + diag[row] + 1; 19478862d2efSBarry Smith idx = a->j + diag[row] + 1; 19482af78befSBarry Smith 19498862d2efSBarry Smith /* see comments for MatMult_Inode() for how this is coded */ 19508862d2efSBarry Smith switch (sizes[i]){ 19518862d2efSBarry Smith case 1: 19528862d2efSBarry Smith 19538862d2efSBarry Smith sum1 = xb[row]; 19548862d2efSBarry Smith for(n = 0; n<sz-1; n+=2) { 19558862d2efSBarry Smith i1 = idx[0]; 19568862d2efSBarry Smith i2 = idx[1]; 19578862d2efSBarry Smith idx += 2; 19588862d2efSBarry Smith tmp0 = x[i1]; 19598862d2efSBarry Smith tmp1 = x[i2]; 19608862d2efSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 19618862d2efSBarry Smith } 19628862d2efSBarry Smith 19638862d2efSBarry Smith if (n == sz-1){ 1964f0d39aaaSBarry Smith tmp0 = x[*idx]; 1965f0d39aaaSBarry Smith sum1 -= *v1*tmp0; 19668862d2efSBarry Smith } 1967f0d39aaaSBarry Smith x[row--] = sum1*(*ibdiag); 1968f0d39aaaSBarry Smith break; 1969f0d39aaaSBarry Smith 1970f0d39aaaSBarry Smith case 2: 1971f0d39aaaSBarry Smith 1972f0d39aaaSBarry Smith sum1 = xb[row]; 1973f0d39aaaSBarry Smith sum2 = xb[row-1]; 1974f0d39aaaSBarry Smith /* note that sum1 is associated with the second of the two rows */ 1975f0d39aaaSBarry Smith v2 = a->a + diag[row-1] + 2; 1976f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 1977f0d39aaaSBarry Smith i1 = idx[0]; 1978f0d39aaaSBarry Smith i2 = idx[1]; 1979f0d39aaaSBarry Smith idx += 2; 1980f0d39aaaSBarry Smith tmp0 = x[i1]; 1981f0d39aaaSBarry Smith tmp1 = x[i2]; 1982f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 1983f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 1984f0d39aaaSBarry Smith } 1985f0d39aaaSBarry Smith 1986f0d39aaaSBarry Smith if (n == sz-1){ 1987f0d39aaaSBarry Smith tmp0 = x[*idx]; 1988f0d39aaaSBarry Smith sum1 -= *v1*tmp0; 1989f0d39aaaSBarry Smith sum2 -= *v2*tmp0; 1990f0d39aaaSBarry Smith } 1991f0d39aaaSBarry Smith x[row--] = sum2*ibdiag[1] + sum1*ibdiag[3]; 1992f0d39aaaSBarry Smith x[row--] = sum2*ibdiag[0] + sum1*ibdiag[2]; 1993f0d39aaaSBarry Smith break; 1994f0d39aaaSBarry Smith case 3: 1995f0d39aaaSBarry Smith 1996f0d39aaaSBarry Smith sum1 = xb[row]; 1997f0d39aaaSBarry Smith sum2 = xb[row-1]; 1998f0d39aaaSBarry Smith sum3 = xb[row-2]; 1999f0d39aaaSBarry Smith v2 = a->a + diag[row-1] + 2; 2000f0d39aaaSBarry Smith v3 = a->a + diag[row-2] + 3; 2001f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 2002f0d39aaaSBarry Smith i1 = idx[0]; 2003f0d39aaaSBarry Smith i2 = idx[1]; 2004f0d39aaaSBarry Smith idx += 2; 2005f0d39aaaSBarry Smith tmp0 = x[i1]; 2006f0d39aaaSBarry Smith tmp1 = x[i2]; 2007f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 2008f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 2009f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 2010f0d39aaaSBarry Smith } 2011f0d39aaaSBarry Smith 2012f0d39aaaSBarry Smith if (n == sz-1){ 2013f0d39aaaSBarry Smith tmp0 = x[*idx]; 2014f0d39aaaSBarry Smith sum1 -= *v1*tmp0; 2015f0d39aaaSBarry Smith sum2 -= *v2*tmp0; 2016f0d39aaaSBarry Smith sum3 -= *v3*tmp0; 2017f0d39aaaSBarry Smith } 2018f0d39aaaSBarry Smith x[row--] = sum3*ibdiag[2] + sum2*ibdiag[5] + sum1*ibdiag[8]; 2019f0d39aaaSBarry Smith x[row--] = sum3*ibdiag[1] + sum2*ibdiag[4] + sum1*ibdiag[7]; 2020f0d39aaaSBarry Smith x[row--] = sum3*ibdiag[0] + sum2*ibdiag[3] + sum1*ibdiag[6]; 2021f0d39aaaSBarry Smith break; 2022f0d39aaaSBarry Smith case 4: 2023f0d39aaaSBarry Smith 2024f0d39aaaSBarry Smith sum1 = xb[row]; 2025f0d39aaaSBarry Smith sum2 = xb[row-1]; 2026f0d39aaaSBarry Smith sum3 = xb[row-2]; 2027f0d39aaaSBarry Smith sum4 = xb[row-3]; 2028f0d39aaaSBarry Smith v2 = a->a + diag[row-1] + 2; 2029f0d39aaaSBarry Smith v3 = a->a + diag[row-2] + 3; 2030f0d39aaaSBarry Smith v4 = a->a + diag[row-3] + 4; 2031f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 2032f0d39aaaSBarry Smith i1 = idx[0]; 2033f0d39aaaSBarry Smith i2 = idx[1]; 2034f0d39aaaSBarry Smith idx += 2; 2035f0d39aaaSBarry Smith tmp0 = x[i1]; 2036f0d39aaaSBarry Smith tmp1 = x[i2]; 2037f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 2038f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 2039f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 2040f0d39aaaSBarry Smith sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2; 2041f0d39aaaSBarry Smith } 2042f0d39aaaSBarry Smith 2043f0d39aaaSBarry Smith if (n == sz-1){ 2044f0d39aaaSBarry Smith tmp0 = x[*idx]; 2045f0d39aaaSBarry Smith sum1 -= *v1*tmp0; 2046f0d39aaaSBarry Smith sum2 -= *v2*tmp0; 2047f0d39aaaSBarry Smith sum3 -= *v3*tmp0; 2048f0d39aaaSBarry Smith sum4 -= *v4*tmp0; 2049f0d39aaaSBarry Smith } 2050f0d39aaaSBarry Smith x[row--] = sum4*ibdiag[3] + sum3*ibdiag[7] + sum2*ibdiag[11] + sum1*ibdiag[15]; 2051f0d39aaaSBarry Smith x[row--] = sum4*ibdiag[2] + sum3*ibdiag[6] + sum2*ibdiag[10] + sum1*ibdiag[14]; 2052f0d39aaaSBarry Smith x[row--] = sum4*ibdiag[1] + sum3*ibdiag[5] + sum2*ibdiag[9] + sum1*ibdiag[13]; 2053f0d39aaaSBarry Smith x[row--] = sum4*ibdiag[0] + sum3*ibdiag[4] + sum2*ibdiag[8] + sum1*ibdiag[12]; 2054f0d39aaaSBarry Smith break; 2055f0d39aaaSBarry Smith case 5: 2056f0d39aaaSBarry Smith 2057f0d39aaaSBarry Smith sum1 = xb[row]; 2058f0d39aaaSBarry Smith sum2 = xb[row-1]; 2059f0d39aaaSBarry Smith sum3 = xb[row-2]; 2060f0d39aaaSBarry Smith sum4 = xb[row-3]; 2061f0d39aaaSBarry Smith sum5 = xb[row-4]; 2062f0d39aaaSBarry Smith v2 = a->a + diag[row-1] + 2; 2063f0d39aaaSBarry Smith v3 = a->a + diag[row-2] + 3; 2064f0d39aaaSBarry Smith v4 = a->a + diag[row-3] + 4; 2065f0d39aaaSBarry Smith v5 = a->a + diag[row-4] + 5; 2066f0d39aaaSBarry Smith for(n = 0; n<sz-1; n+=2) { 2067f0d39aaaSBarry Smith i1 = idx[0]; 2068f0d39aaaSBarry Smith i2 = idx[1]; 2069f0d39aaaSBarry Smith idx += 2; 2070f0d39aaaSBarry Smith tmp0 = x[i1]; 2071f0d39aaaSBarry Smith tmp1 = x[i2]; 2072f0d39aaaSBarry Smith sum1 -= v1[0] * tmp0 + v1[1] * tmp1; v1 += 2; 2073f0d39aaaSBarry Smith sum2 -= v2[0] * tmp0 + v2[1] * tmp1; v2 += 2; 2074f0d39aaaSBarry Smith sum3 -= v3[0] * tmp0 + v3[1] * tmp1; v3 += 2; 2075f0d39aaaSBarry Smith sum4 -= v4[0] * tmp0 + v4[1] * tmp1; v4 += 2; 2076f0d39aaaSBarry Smith sum5 -= v5[0] * tmp0 + v5[1] * tmp1; v5 += 2; 2077f0d39aaaSBarry Smith } 2078f0d39aaaSBarry Smith 2079f0d39aaaSBarry Smith if (n == sz-1){ 2080f0d39aaaSBarry Smith tmp0 = x[*idx]; 2081f0d39aaaSBarry Smith sum1 -= *v1*tmp0; 2082f0d39aaaSBarry Smith sum2 -= *v2*tmp0; 2083f0d39aaaSBarry Smith sum3 -= *v3*tmp0; 2084f0d39aaaSBarry Smith sum4 -= *v4*tmp0; 2085f0d39aaaSBarry Smith sum5 -= *v5*tmp0; 2086f0d39aaaSBarry Smith } 2087f0d39aaaSBarry Smith x[row--] = sum5*ibdiag[4] + sum4*ibdiag[9] + sum3*ibdiag[14] + sum2*ibdiag[19] + sum1*ibdiag[24]; 2088f0d39aaaSBarry Smith x[row--] = sum5*ibdiag[3] + sum4*ibdiag[8] + sum3*ibdiag[13] + sum2*ibdiag[18] + sum1*ibdiag[23]; 2089f0d39aaaSBarry Smith x[row--] = sum5*ibdiag[2] + sum4*ibdiag[7] + sum3*ibdiag[12] + sum2*ibdiag[17] + sum1*ibdiag[22]; 2090f0d39aaaSBarry Smith x[row--] = sum5*ibdiag[1] + sum4*ibdiag[6] + sum3*ibdiag[11] + sum2*ibdiag[16] + sum1*ibdiag[21]; 2091f0d39aaaSBarry Smith x[row--] = sum5*ibdiag[0] + sum4*ibdiag[5] + sum3*ibdiag[10] + sum2*ibdiag[15] + sum1*ibdiag[20]; 20928862d2efSBarry Smith break; 20938862d2efSBarry Smith default: 20948862d2efSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Inode size %D not supported",sizes[i]); 20958862d2efSBarry Smith } 20962af78befSBarry Smith } 20972af78befSBarry Smith 20982af78befSBarry Smith ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20992af78befSBarry Smith } 21002af78befSBarry Smith its--; 21012af78befSBarry Smith } 21028862d2efSBarry Smith if (its) SETERRQ(PETSC_ERR_SUP,"Currently no support for multiply SOR sweeps using inode version of AIJ matrix format;\n run with the option -mat_no_inode"); 21032af78befSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 21042af78befSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);} 21052af78befSBarry Smith PetscFunctionReturn(0); 21062af78befSBarry Smith } 21072af78befSBarry Smith 21082af78befSBarry Smith 21094c1414c8SBarry Smith /* 21104c1414c8SBarry Smith samestructure indicates that the matrix has not changed its nonzero structure so we 21114c1414c8SBarry Smith do not need to recompute the inodes 21124c1414c8SBarry Smith */ 21134c1414c8SBarry Smith #undef __FUNCT__ 21144c1414c8SBarry Smith #define __FUNCT__ "Mat_CheckInode" 21154c1414c8SBarry Smith PetscErrorCode Mat_CheckInode(Mat A,PetscTruth samestructure) 21164c1414c8SBarry Smith { 21174c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 21184c1414c8SBarry Smith PetscErrorCode ierr; 21194c1414c8SBarry Smith PetscInt i,j,m,nzx,nzy,*idx,*idy,*ns,*ii,node_count,blk_size; 2120d74fe821SBarry Smith PetscTruth flag; 21214c1414c8SBarry Smith 21224c1414c8SBarry Smith PetscFunctionBegin; 2123d74fe821SBarry Smith if (!a->inode.use) PetscFunctionReturn(0); 21244c1414c8SBarry Smith if (a->inode.checked && samestructure) PetscFunctionReturn(0); 21254c1414c8SBarry Smith 21264c1414c8SBarry Smith 2127d0f46423SBarry Smith m = A->rmap->n; 21284c1414c8SBarry Smith if (a->inode.size) {ns = a->inode.size;} 21294c1414c8SBarry Smith else {ierr = PetscMalloc((m+1)*sizeof(PetscInt),&ns);CHKERRQ(ierr);} 21304c1414c8SBarry Smith 21314c1414c8SBarry Smith i = 0; 21324c1414c8SBarry Smith node_count = 0; 21334c1414c8SBarry Smith idx = a->j; 21344c1414c8SBarry Smith ii = a->i; 21354c1414c8SBarry Smith while (i < m){ /* For each row */ 21364c1414c8SBarry Smith nzx = ii[i+1] - ii[i]; /* Number of nonzeros */ 21374c1414c8SBarry Smith /* Limits the number of elements in a node to 'a->inode.limit' */ 21384c1414c8SBarry Smith for (j=i+1,idy=idx,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) { 21394c1414c8SBarry Smith nzy = ii[j+1] - ii[j]; /* Same number of nonzeros */ 21404c1414c8SBarry Smith if (nzy != nzx) break; 21414c1414c8SBarry Smith idy += nzx; /* Same nonzero pattern */ 21424c1414c8SBarry Smith ierr = PetscMemcmp(idx,idy,nzx*sizeof(PetscInt),&flag);CHKERRQ(ierr); 21434c1414c8SBarry Smith if (!flag) break; 21444c1414c8SBarry Smith } 21454c1414c8SBarry Smith ns[node_count++] = blk_size; 21464c1414c8SBarry Smith idx += blk_size*nzx; 21474c1414c8SBarry Smith i = j; 21484c1414c8SBarry Smith } 21494c1414c8SBarry Smith /* If not enough inodes found,, do not use inode version of the routines */ 2150f0d39aaaSBarry Smith if (!a->inode.size && m && node_count > .9*m) { 21514c1414c8SBarry Smith ierr = PetscFree(ns);CHKERRQ(ierr); 21524c1414c8SBarry Smith a->inode.node_count = 0; 21534c1414c8SBarry Smith a->inode.size = PETSC_NULL; 21544c1414c8SBarry Smith a->inode.use = PETSC_FALSE; 21554c1414c8SBarry Smith ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr); 21564c1414c8SBarry Smith } else { 21574c1414c8SBarry Smith A->ops->mult = MatMult_Inode; 21582af78befSBarry Smith A->ops->relax = MatRelax_Inode; 21594c1414c8SBarry Smith A->ops->multadd = MatMultAdd_Inode; 21604c1414c8SBarry Smith A->ops->getrowij = MatGetRowIJ_Inode; 21614c1414c8SBarry Smith A->ops->restorerowij = MatRestoreRowIJ_Inode; 21624c1414c8SBarry Smith A->ops->getcolumnij = MatGetColumnIJ_Inode; 21634c1414c8SBarry Smith A->ops->restorecolumnij = MatRestoreColumnIJ_Inode; 21644c1414c8SBarry Smith A->ops->coloringpatch = MatColoringPatch_Inode; 21654c1414c8SBarry Smith a->inode.node_count = node_count; 21664c1414c8SBarry Smith a->inode.size = ns; 21674c1414c8SBarry Smith ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr); 21684c1414c8SBarry Smith } 21694c1414c8SBarry Smith PetscFunctionReturn(0); 21704c1414c8SBarry Smith } 21714c1414c8SBarry Smith 21724c1414c8SBarry Smith /* 21734c1414c8SBarry Smith This is really ugly. if inodes are used this replaces the 21744c1414c8SBarry Smith permutations with ones that correspond to rows/cols of the matrix 21754c1414c8SBarry Smith rather then inode blocks 21764c1414c8SBarry Smith */ 21774c1414c8SBarry Smith #undef __FUNCT__ 21784c1414c8SBarry Smith #define __FUNCT__ "MatInodeAdjustForInodes" 21794c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat A,IS *rperm,IS *cperm) 21804c1414c8SBarry Smith { 21814c1414c8SBarry Smith PetscErrorCode ierr,(*f)(Mat,IS*,IS*); 21824c1414c8SBarry Smith 21834c1414c8SBarry Smith PetscFunctionBegin; 21844c1414c8SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatInodeAdjustForInodes_C",(void (**)(void))&f);CHKERRQ(ierr); 21854c1414c8SBarry Smith if (f) { 21864c1414c8SBarry Smith ierr = (*f)(A,rperm,cperm);CHKERRQ(ierr); 21874c1414c8SBarry Smith } 21884c1414c8SBarry Smith PetscFunctionReturn(0); 21894c1414c8SBarry Smith } 21904c1414c8SBarry Smith 21914c1414c8SBarry Smith EXTERN_C_BEGIN 21924c1414c8SBarry Smith #undef __FUNCT__ 21934c1414c8SBarry Smith #define __FUNCT__ "MatAdjustForInodes_Inode" 21944c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes_Inode(Mat A,IS *rperm,IS *cperm) 21954c1414c8SBarry Smith { 21964c1414c8SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 21974c1414c8SBarry Smith PetscErrorCode ierr; 2198d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,i,j,*ridx,*cidx,nslim_row = a->inode.node_count; 21994c1414c8SBarry Smith PetscInt row,col,*permr,*permc,*ns_row = a->inode.size,*tns,start_val,end_val,indx; 22004c1414c8SBarry Smith PetscInt nslim_col,*ns_col; 22014c1414c8SBarry Smith IS ris = *rperm,cis = *cperm; 22024c1414c8SBarry Smith 22034c1414c8SBarry Smith PetscFunctionBegin; 22044c1414c8SBarry Smith if (!a->inode.size) PetscFunctionReturn(0); /* no inodes so return */ 22054c1414c8SBarry Smith if (a->inode.node_count == m) PetscFunctionReturn(0); /* all inodes are of size 1 */ 22064c1414c8SBarry Smith 22074c1414c8SBarry Smith ierr = Mat_CreateColInode(A,&nslim_col,&ns_col);CHKERRQ(ierr); 22084c1414c8SBarry Smith ierr = PetscMalloc((((nslim_row>nslim_col)?nslim_row:nslim_col)+1)*sizeof(PetscInt),&tns);CHKERRQ(ierr); 22094c1414c8SBarry Smith ierr = PetscMalloc((m+n+1)*sizeof(PetscInt),&permr);CHKERRQ(ierr); 22104c1414c8SBarry Smith permc = permr + m; 22114c1414c8SBarry Smith 22124c1414c8SBarry Smith ierr = ISGetIndices(ris,&ridx);CHKERRQ(ierr); 22134c1414c8SBarry Smith ierr = ISGetIndices(cis,&cidx);CHKERRQ(ierr); 22144c1414c8SBarry Smith 22154c1414c8SBarry Smith /* Form the inode structure for the rows of permuted matric using inv perm*/ 22164c1414c8SBarry Smith for (i=0,tns[0]=0; i<nslim_row; ++i) tns[i+1] = tns[i] + ns_row[i]; 22174c1414c8SBarry Smith 22184c1414c8SBarry Smith /* Construct the permutations for rows*/ 22194c1414c8SBarry Smith for (i=0,row = 0; i<nslim_row; ++i){ 22204c1414c8SBarry Smith indx = ridx[i]; 22214c1414c8SBarry Smith start_val = tns[indx]; 22224c1414c8SBarry Smith end_val = tns[indx + 1]; 22234c1414c8SBarry Smith for (j=start_val; j<end_val; ++j,++row) permr[row]= j; 22244c1414c8SBarry Smith } 22254c1414c8SBarry Smith 22264c1414c8SBarry Smith /* Form the inode structure for the columns of permuted matrix using inv perm*/ 22274c1414c8SBarry Smith for (i=0,tns[0]=0; i<nslim_col; ++i) tns[i+1] = tns[i] + ns_col[i]; 22284c1414c8SBarry Smith 22294c1414c8SBarry Smith /* Construct permutations for columns */ 22304c1414c8SBarry Smith for (i=0,col=0; i<nslim_col; ++i){ 22314c1414c8SBarry Smith indx = cidx[i]; 22324c1414c8SBarry Smith start_val = tns[indx]; 22334c1414c8SBarry Smith end_val = tns[indx + 1]; 22344c1414c8SBarry Smith for (j = start_val; j<end_val; ++j,++col) permc[col]= j; 22354c1414c8SBarry Smith } 22364c1414c8SBarry Smith 22374c1414c8SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,permr,rperm);CHKERRQ(ierr); 22384c1414c8SBarry Smith ISSetPermutation(*rperm); 22394c1414c8SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,permc,cperm);CHKERRQ(ierr); 22404c1414c8SBarry Smith ISSetPermutation(*cperm); 22414c1414c8SBarry Smith 22424c1414c8SBarry Smith ierr = ISRestoreIndices(ris,&ridx);CHKERRQ(ierr); 22434c1414c8SBarry Smith ierr = ISRestoreIndices(cis,&cidx);CHKERRQ(ierr); 22444c1414c8SBarry Smith 22454c1414c8SBarry Smith ierr = PetscFree(ns_col);CHKERRQ(ierr); 22464c1414c8SBarry Smith ierr = PetscFree(permr);CHKERRQ(ierr); 22474c1414c8SBarry Smith ierr = ISDestroy(cis);CHKERRQ(ierr); 22484c1414c8SBarry Smith ierr = ISDestroy(ris);CHKERRQ(ierr); 22494c1414c8SBarry Smith ierr = PetscFree(tns);CHKERRQ(ierr); 22504c1414c8SBarry Smith PetscFunctionReturn(0); 22514c1414c8SBarry Smith } 22524c1414c8SBarry Smith EXTERN_C_END 22534c1414c8SBarry Smith 22544c1414c8SBarry Smith #undef __FUNCT__ 22554c1414c8SBarry Smith #define __FUNCT__ "MatInodeGetInodeSizes" 22564c1414c8SBarry Smith /*@C 22574c1414c8SBarry Smith MatInodeGetInodeSizes - Returns the inode information of the Inode matrix. 22584c1414c8SBarry Smith 22594c1414c8SBarry Smith Collective on Mat 22604c1414c8SBarry Smith 22614c1414c8SBarry Smith Input Parameter: 22624c1414c8SBarry Smith . A - the Inode matrix or matrix derived from the Inode class -- e.g., SeqAIJ 22634c1414c8SBarry Smith 22644c1414c8SBarry Smith Output Parameter: 22654c1414c8SBarry Smith + node_count - no of inodes present in the matrix. 22664c1414c8SBarry Smith . sizes - an array of size node_count,with sizes of each inode. 22674c1414c8SBarry Smith - limit - the max size used to generate the inodes. 22684c1414c8SBarry Smith 22694c1414c8SBarry Smith Level: advanced 22704c1414c8SBarry Smith 22714c1414c8SBarry Smith Notes: This routine returns some internal storage information 22724c1414c8SBarry Smith of the matrix, it is intended to be used by advanced users. 22734c1414c8SBarry Smith It should be called after the matrix is assembled. 22744c1414c8SBarry Smith The contents of the sizes[] array should not be changed. 22754c1414c8SBarry Smith PETSC_NULL may be passed for information not requested. 22764c1414c8SBarry Smith 22774c1414c8SBarry Smith .keywords: matrix, seqaij, get, inode 22784c1414c8SBarry Smith 22794c1414c8SBarry Smith .seealso: MatGetInfo() 22804c1414c8SBarry Smith @*/ 22814c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat A,PetscInt *node_count,PetscInt *sizes[],PetscInt *limit) 22824c1414c8SBarry Smith { 22834c1414c8SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*[],PetscInt*); 22844c1414c8SBarry Smith 22854c1414c8SBarry Smith PetscFunctionBegin; 22864c1414c8SBarry Smith if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 22874c1414c8SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatInodeGetInodeSizes_C",(void (**)(void))&f);CHKERRQ(ierr); 22884c1414c8SBarry Smith if (f) { 22894c1414c8SBarry Smith ierr = (*f)(A,node_count,sizes,limit);CHKERRQ(ierr); 22904c1414c8SBarry Smith } 22914c1414c8SBarry Smith PetscFunctionReturn(0); 22924c1414c8SBarry Smith } 22934c1414c8SBarry Smith 22944c1414c8SBarry Smith EXTERN_C_BEGIN 22954c1414c8SBarry Smith #undef __FUNCT__ 22964c1414c8SBarry Smith #define __FUNCT__ "MatInodeGetInodeSizes_Inode" 22974c1414c8SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes_Inode(Mat A,PetscInt *node_count,PetscInt *sizes[],PetscInt *limit) 22984c1414c8SBarry Smith { 22994c1414c8SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 23004c1414c8SBarry Smith 23014c1414c8SBarry Smith PetscFunctionBegin; 23024c1414c8SBarry Smith if (node_count) *node_count = a->inode.node_count; 23034c1414c8SBarry Smith if (sizes) *sizes = a->inode.size; 23044c1414c8SBarry Smith if (limit) *limit = a->inode.limit; 23054c1414c8SBarry Smith PetscFunctionReturn(0); 23064c1414c8SBarry Smith } 23074c1414c8SBarry Smith EXTERN_C_END 2308