1b6490206SBarry Smith 22593348eSBarry Smith #ifndef lint 3*bcd2baecSBarry Smith static char vcid[] = "$Id: baij.c,v 1.9 1996/03/04 05:16:18 bsmith Exp bsmith $"; 42593348eSBarry Smith #endif 52593348eSBarry Smith 62593348eSBarry Smith /* 7b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 82593348eSBarry Smith matrix storage format. 92593348eSBarry Smith */ 10b6490206SBarry Smith #include "baij.h" 112593348eSBarry Smith #include "vec/vecimpl.h" 122593348eSBarry Smith #include "inline/spops.h" 132593348eSBarry Smith #include "petsc.h" 142593348eSBarry Smith 15*bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 162593348eSBarry Smith 17b6490206SBarry Smith static int MatGetReordering_SeqBAIJ(Mat A,MatOrdering type,IS *rperm,IS *cperm) 182593348eSBarry Smith { 19b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 20*bcd2baecSBarry Smith int ierr, *ia, *ja,n = a->mbs,*idx,i,ishift,oshift; 212593348eSBarry Smith 222593348eSBarry Smith /* 232593348eSBarry Smith this is tacky: In the future when we have written special factorization 242593348eSBarry Smith and solve routines for the identity permutation we should use a 252593348eSBarry Smith stride index set instead of the general one. 262593348eSBarry Smith */ 272593348eSBarry Smith if (type == ORDER_NATURAL) { 282593348eSBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 292593348eSBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 302593348eSBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 312593348eSBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 322593348eSBarry Smith PetscFree(idx); 332593348eSBarry Smith ISSetPermutation(*rperm); 342593348eSBarry Smith ISSetPermutation(*cperm); 352593348eSBarry Smith ISSetIdentity(*rperm); 362593348eSBarry Smith ISSetIdentity(*cperm); 372593348eSBarry Smith return 0; 382593348eSBarry Smith } 392593348eSBarry Smith 40*bcd2baecSBarry Smith MatReorderingRegisterAll(); 41*bcd2baecSBarry Smith ishift = a->indexshift; 42*bcd2baecSBarry Smith oshift = -MatReorderingIndexShift[(int)type]; 43*bcd2baecSBarry Smith if (MatReorderingRequiresSymmetric[(int)type]) { 44*bcd2baecSBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja); 45*bcd2baecSBarry Smith CHKERRQ(ierr); 46*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 472593348eSBarry Smith PetscFree(ia); PetscFree(ja); 48*bcd2baecSBarry Smith } else { 49*bcd2baecSBarry Smith if (ishift == oshift) { 50*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 51*bcd2baecSBarry Smith } 52*bcd2baecSBarry Smith else if (ishift == -1) { 53*bcd2baecSBarry Smith /* temporarily subtract 1 from i and j indices */ 54*bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 55*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 56*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 57*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 58*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 59*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 60*bcd2baecSBarry Smith } else { 61*bcd2baecSBarry Smith /* temporarily add 1 to i and j indices */ 62*bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 63*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 64*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 65*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 66*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 67*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 68*bcd2baecSBarry Smith } 69*bcd2baecSBarry Smith } 702593348eSBarry Smith return 0; 712593348eSBarry Smith } 722593348eSBarry Smith 73de6a44a3SBarry Smith /* 74de6a44a3SBarry Smith Adds diagonal pointers to sparse matrix structure. 75de6a44a3SBarry Smith */ 76de6a44a3SBarry Smith 77de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 78de6a44a3SBarry Smith { 79de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 807fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 81de6a44a3SBarry Smith 82de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 83de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 847fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 85de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 86de6a44a3SBarry Smith if (a->j[j] == i) { 87de6a44a3SBarry Smith diag[i] = j; 88de6a44a3SBarry Smith break; 89de6a44a3SBarry Smith } 90de6a44a3SBarry Smith } 91de6a44a3SBarry Smith } 92de6a44a3SBarry Smith a->diag = diag; 93de6a44a3SBarry Smith return 0; 94de6a44a3SBarry Smith } 952593348eSBarry Smith 962593348eSBarry Smith #include "draw.h" 972593348eSBarry Smith #include "pinclude/pviewer.h" 982593348eSBarry Smith #include "sysio.h" 992593348eSBarry Smith 100b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 1012593348eSBarry Smith { 102b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 103b6490206SBarry Smith int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l; 104b6490206SBarry Smith Scalar *aa; 1052593348eSBarry Smith 1062593348eSBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1072593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 1082593348eSBarry Smith col_lens[0] = MAT_COOKIE; 1092593348eSBarry Smith col_lens[1] = a->m; 1102593348eSBarry Smith col_lens[2] = a->n; 111b6490206SBarry Smith col_lens[3] = a->nz*bs*bs; 1122593348eSBarry Smith 1132593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 114b6490206SBarry Smith count = 0; 115b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 116b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 117b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 118b6490206SBarry Smith } 1192593348eSBarry Smith } 1202593348eSBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 1212593348eSBarry Smith PetscFree(col_lens); 1222593348eSBarry Smith 1232593348eSBarry Smith /* store column indices (zero start index) */ 124b6490206SBarry Smith jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj); 125b6490206SBarry Smith count = 0; 126b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 127b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 128b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 129b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 130b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 1312593348eSBarry Smith } 1322593348eSBarry Smith } 133b6490206SBarry Smith } 134b6490206SBarry Smith } 135b6490206SBarry Smith ierr = SYWrite(fd,jj,bs*bs*a->nz,SYINT,0); CHKERRQ(ierr); 136b6490206SBarry Smith PetscFree(jj); 1372593348eSBarry Smith 1382593348eSBarry Smith /* store nonzero values */ 139b6490206SBarry Smith aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa); 140b6490206SBarry Smith count = 0; 141b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 142b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 143b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 144b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 145b6490206SBarry Smith aa[count++] = a->a[bs*bs*k + l*bs + j]; 146b6490206SBarry Smith } 147b6490206SBarry Smith } 148b6490206SBarry Smith } 149b6490206SBarry Smith } 150b6490206SBarry Smith ierr = SYWrite(fd,aa,bs*bs*a->nz,SYSCALAR,0); CHKERRQ(ierr); 151b6490206SBarry Smith PetscFree(aa); 1522593348eSBarry Smith return 0; 1532593348eSBarry Smith } 1542593348eSBarry Smith 155b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 1562593348eSBarry Smith { 157b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 158de6a44a3SBarry Smith int ierr, i,j,format,bs = a->bs,k,l; 1592593348eSBarry Smith FILE *fd; 1602593348eSBarry Smith char *outputname; 1612593348eSBarry Smith 1622593348eSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 1632593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 1642593348eSBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 1652593348eSBarry Smith if (format == FILE_FORMAT_INFO) { 1662593348eSBarry Smith /* no need to print additional information */ ; 1672593348eSBarry Smith } 1682593348eSBarry Smith else if (format == FILE_FORMAT_MATLAB) { 169b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported"); 1702593348eSBarry Smith } 1712593348eSBarry Smith else { 172b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 173b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 174b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 175b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 176b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 17788685aaeSLois Curfman McInnes #if defined(PETSC_COMPLEX) 17888685aaeSLois Curfman McInnes if (imag(a->a[bs*bs*k + l*bs + j]) != 0.0) { 17988685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 18088685aaeSLois Curfman McInnes real(a->a[bs*bs*k + l*bs + j]),imag(a->a[bs*bs*k + l*bs + j])); 18188685aaeSLois Curfman McInnes } 18288685aaeSLois Curfman McInnes else { 18388685aaeSLois Curfman McInnes fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs*bs*k + l*bs + j])); 18488685aaeSLois Curfman McInnes } 18588685aaeSLois Curfman McInnes #else 186de6a44a3SBarry Smith fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs*bs*k + l*bs + j]); 18788685aaeSLois Curfman McInnes #endif 1882593348eSBarry Smith } 1892593348eSBarry Smith } 1902593348eSBarry Smith fprintf(fd,"\n"); 1912593348eSBarry Smith } 1922593348eSBarry Smith } 193b6490206SBarry Smith } 1942593348eSBarry Smith fflush(fd); 1952593348eSBarry Smith return 0; 1962593348eSBarry Smith } 1972593348eSBarry Smith 198b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer) 1992593348eSBarry Smith { 2002593348eSBarry Smith Mat A = (Mat) obj; 2012593348eSBarry Smith PetscObject vobj = (PetscObject) viewer; 2022593348eSBarry Smith 2032593348eSBarry Smith if (!viewer) { 2042593348eSBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 2052593348eSBarry Smith } 2062593348eSBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 2072593348eSBarry Smith if (vobj->type == MATLAB_VIEWER) { 208b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported"); 2092593348eSBarry Smith } 2102593348eSBarry Smith else if (vobj->type == ASCII_FILE_VIEWER||vobj->type == ASCII_FILES_VIEWER){ 211b6490206SBarry Smith return MatView_SeqBAIJ_ASCII(A,viewer); 2122593348eSBarry Smith } 2132593348eSBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 214b6490206SBarry Smith return MatView_SeqBAIJ_Binary(A,viewer); 2152593348eSBarry Smith } 2162593348eSBarry Smith } 2172593348eSBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 2182593348eSBarry Smith if (vobj->type == NULLWINDOW) return 0; 219b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported"); 2202593348eSBarry Smith } 2212593348eSBarry Smith return 0; 2222593348eSBarry Smith } 223b6490206SBarry Smith 224b6490206SBarry Smith 225b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A) 2262593348eSBarry Smith { 227b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 228de6a44a3SBarry Smith PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar)); 2292593348eSBarry Smith return 0; 2302593348eSBarry Smith } 2312593348eSBarry Smith 232b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj) 2332593348eSBarry Smith { 2342593348eSBarry Smith Mat A = (Mat) obj; 235b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2362593348eSBarry Smith 2372593348eSBarry Smith #if defined(PETSC_LOG) 2382593348eSBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 2392593348eSBarry Smith #endif 2402593348eSBarry Smith PetscFree(a->a); 2412593348eSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 2422593348eSBarry Smith if (a->diag) PetscFree(a->diag); 2432593348eSBarry Smith if (a->ilen) PetscFree(a->ilen); 2442593348eSBarry Smith if (a->imax) PetscFree(a->imax); 2452593348eSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 246de6a44a3SBarry Smith if (a->mult_work) PetscFree(a->mult_work); 2472593348eSBarry Smith PetscFree(a); 2482593348eSBarry Smith PLogObjectDestroy(A); 2492593348eSBarry Smith PetscHeaderDestroy(A); 2502593348eSBarry Smith return 0; 2512593348eSBarry Smith } 2522593348eSBarry Smith 253b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op) 2542593348eSBarry Smith { 255b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2562593348eSBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 2572593348eSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 2582593348eSBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 2592593348eSBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 2602593348eSBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 2612593348eSBarry Smith else if (op == ROWS_SORTED || 2622593348eSBarry Smith op == SYMMETRIC_MATRIX || 2632593348eSBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 2642593348eSBarry Smith op == YES_NEW_DIAGONALS) 265b6490206SBarry Smith PLogInfo((PetscObject)A,"Info:MatSetOption_SeqBAIJ:Option ignored\n"); 2662593348eSBarry Smith else if (op == NO_NEW_DIAGONALS) 267b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");} 2682593348eSBarry Smith else 269b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");} 2702593348eSBarry Smith return 0; 2712593348eSBarry Smith } 2722593348eSBarry Smith 2732593348eSBarry Smith 2742593348eSBarry Smith /* -------------------------------------------------------*/ 2752593348eSBarry Smith /* Should check that shapes of vectors and matrices match */ 2762593348eSBarry Smith /* -------------------------------------------------------*/ 277b6490206SBarry Smith #include "pinclude/plapack.h" 278b6490206SBarry Smith 279b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy) 2802593348eSBarry Smith { 281b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 282b6490206SBarry Smith Scalar *xg,*yg; 283b6490206SBarry Smith register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5; 284b6490206SBarry Smith register Scalar x1,x2,x3,x4,x5; 285b6490206SBarry Smith int mbs = a->mbs, m = a->m, i, *idx,*ii; 286b6490206SBarry Smith int bs = a->bs,j,n,bs2 = bs*bs; 2872593348eSBarry Smith 288b6490206SBarry Smith VecGetArray(xx,&xg); x = xg; VecGetArray(yy,&yg); y = yg; 289b6490206SBarry Smith PetscMemzero(y,m*sizeof(Scalar)); 290b6490206SBarry Smith x = x; 2912593348eSBarry Smith idx = a->j; 2922593348eSBarry Smith v = a->a; 2932593348eSBarry Smith ii = a->i; 294b6490206SBarry Smith 295b6490206SBarry Smith switch (bs) { 296b6490206SBarry Smith case 1: 2972593348eSBarry Smith for ( i=0; i<m; i++ ) { 2982593348eSBarry Smith n = ii[1] - ii[0]; ii++; 2992593348eSBarry Smith sum = 0.0; 3002593348eSBarry Smith while (n--) sum += *v++ * x[*idx++]; 3012593348eSBarry Smith y[i] = sum; 3022593348eSBarry Smith } 3032593348eSBarry Smith break; 304b6490206SBarry Smith case 2: 305b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 306de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 307b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; 308b6490206SBarry Smith for ( j=0; j<n; j++ ) { 309b6490206SBarry Smith xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 310b6490206SBarry Smith sum1 += v[0]*x1 + v[2]*x2; 311b6490206SBarry Smith sum2 += v[1]*x1 + v[3]*x2; 312b6490206SBarry Smith v += 4; 313b6490206SBarry Smith } 314b6490206SBarry Smith y[0] += sum1; y[1] += sum2; 315b6490206SBarry Smith y += 2; 316b6490206SBarry Smith } 317b6490206SBarry Smith break; 318b6490206SBarry Smith case 3: 319b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 320de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 321b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 322b6490206SBarry Smith for ( j=0; j<n; j++ ) { 323b6490206SBarry Smith xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 324b6490206SBarry Smith sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 325b6490206SBarry Smith sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 326b6490206SBarry Smith sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 327b6490206SBarry Smith v += 9; 328b6490206SBarry Smith } 329b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; 330b6490206SBarry Smith y += 3; 331b6490206SBarry Smith } 332b6490206SBarry Smith break; 333b6490206SBarry Smith case 4: 334b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 335de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 336b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 337b6490206SBarry Smith for ( j=0; j<n; j++ ) { 338b6490206SBarry Smith xb = x + 4*(*idx++); 339b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 340b6490206SBarry Smith sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 341b6490206SBarry Smith sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 342b6490206SBarry Smith sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 343b6490206SBarry Smith sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 344b6490206SBarry Smith v += 16; 345b6490206SBarry Smith } 346b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; 347b6490206SBarry Smith y += 4; 348b6490206SBarry Smith } 349b6490206SBarry Smith break; 350b6490206SBarry Smith case 5: 351b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 352de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 353b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 354b6490206SBarry Smith for ( j=0; j<n; j++ ) { 355b6490206SBarry Smith xb = x + 5*(*idx++); 356b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 357b6490206SBarry Smith sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 358b6490206SBarry Smith sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 359b6490206SBarry Smith sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 360b6490206SBarry Smith sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 361b6490206SBarry Smith sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 362b6490206SBarry Smith v += 25; 363b6490206SBarry Smith } 364b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5; 365b6490206SBarry Smith y += 5; 366b6490206SBarry Smith } 367b6490206SBarry Smith break; 368b6490206SBarry Smith /* block sizes larger then 5 by 5 are handled by BLAS */ 369b6490206SBarry Smith default: { 370de6a44a3SBarry Smith int _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt; 371de6a44a3SBarry Smith if (!a->mult_work) { 372de6a44a3SBarry Smith a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar)); 373de6a44a3SBarry Smith CHKPTRQ(a->mult_work); 374de6a44a3SBarry Smith } 375de6a44a3SBarry Smith work = a->mult_work; 376b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 377de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 378de6a44a3SBarry Smith ncols = n*bs; 379de6a44a3SBarry Smith workt = work; 380b6490206SBarry Smith for ( j=0; j<n; j++ ) { 381b6490206SBarry Smith xb = x + bs*(*idx++); 382de6a44a3SBarry Smith for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 383de6a44a3SBarry Smith workt += bs; 384b6490206SBarry Smith } 385de6a44a3SBarry Smith LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One); 386de6a44a3SBarry Smith v += n*bs2; 387b6490206SBarry Smith y += bs; 3882593348eSBarry Smith } 3892593348eSBarry Smith } 3902593348eSBarry Smith } 391b6490206SBarry Smith PLogFlops(2*bs2*a->nz - m); 3922593348eSBarry Smith return 0; 3932593348eSBarry Smith } 3942593348eSBarry Smith 395de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem) 3962593348eSBarry Smith { 397b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 398*bcd2baecSBarry Smith if (nz) *nz = a->bs*a->bs*a->nz; 399*bcd2baecSBarry Smith if (nza) *nza = a->maxnz; 400*bcd2baecSBarry Smith if (mem) *mem = (int)A->mem; 4012593348eSBarry Smith return 0; 4022593348eSBarry Smith } 4032593348eSBarry Smith 404b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 405b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 406b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec); 4077fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 4087fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 4097fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 4107fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 4117fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 4127fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 4132593348eSBarry Smith 414b6490206SBarry Smith static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 4152593348eSBarry Smith { 416b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4172593348eSBarry Smith *m = a->m; *n = a->n; 4182593348eSBarry Smith return 0; 4192593348eSBarry Smith } 4202593348eSBarry Smith 421b6490206SBarry Smith static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 4222593348eSBarry Smith { 423b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4242593348eSBarry Smith *m = 0; *n = a->m; 4252593348eSBarry Smith return 0; 4262593348eSBarry Smith } 427b6490206SBarry Smith 428b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 4292593348eSBarry Smith { 430b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4312593348eSBarry Smith Scalar *v = a->a; 4322593348eSBarry Smith double sum = 0.0; 433b6490206SBarry Smith int i; 4342593348eSBarry Smith 4352593348eSBarry Smith if (type == NORM_FROBENIUS) { 4362593348eSBarry Smith for (i=0; i<a->nz; i++ ) { 4372593348eSBarry Smith #if defined(PETSC_COMPLEX) 4382593348eSBarry Smith sum += real(conj(*v)*(*v)); v++; 4392593348eSBarry Smith #else 4402593348eSBarry Smith sum += (*v)*(*v); v++; 4412593348eSBarry Smith #endif 4422593348eSBarry Smith } 4432593348eSBarry Smith *norm = sqrt(sum); 4442593348eSBarry Smith } 4452593348eSBarry Smith else { 446b6490206SBarry Smith SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet"); 4472593348eSBarry Smith } 4482593348eSBarry Smith return 0; 4492593348eSBarry Smith } 4502593348eSBarry Smith 4512593348eSBarry Smith /* 4522593348eSBarry Smith note: This can only work for identity for row and col. It would 4532593348eSBarry Smith be good to check this and otherwise generate an error. 4542593348eSBarry Smith */ 455b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 4562593348eSBarry Smith { 457b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 4582593348eSBarry Smith Mat outA; 459de6a44a3SBarry Smith int ierr; 4602593348eSBarry Smith 461b6490206SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported"); 4622593348eSBarry Smith 4632593348eSBarry Smith outA = inA; 4642593348eSBarry Smith inA->factor = FACTOR_LU; 4652593348eSBarry Smith a->row = row; 4662593348eSBarry Smith a->col = col; 4672593348eSBarry Smith 468de6a44a3SBarry Smith a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 4692593348eSBarry Smith 4702593348eSBarry Smith if (!a->diag) { 471b6490206SBarry Smith ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 4722593348eSBarry Smith } 4737fc0212eSBarry Smith ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 4742593348eSBarry Smith return 0; 4752593348eSBarry Smith } 4762593348eSBarry Smith 477b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 4782593348eSBarry Smith { 479b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 480b6490206SBarry Smith int one = 1, totalnz = a->bs*a->bs*a->nz; 481b6490206SBarry Smith BLscal_( &totalnz, alpha, a->a, &one ); 482b6490206SBarry Smith PLogFlops(totalnz); 4832593348eSBarry Smith return 0; 4842593348eSBarry Smith } 4852593348eSBarry Smith 486b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A) 4872593348eSBarry Smith { 4882593348eSBarry Smith static int called = 0; 4892593348eSBarry Smith 4902593348eSBarry Smith if (called) return 0; else called = 1; 4912593348eSBarry Smith return 0; 4922593348eSBarry Smith } 4932593348eSBarry Smith /* -------------------------------------------------------------------*/ 494b6490206SBarry Smith static struct _MatOps MatOps = {0, 495b6490206SBarry Smith 0,0, 496b6490206SBarry Smith MatMult_SeqBAIJ,0, 497b6490206SBarry Smith 0,0, 498de6a44a3SBarry Smith MatSolve_SeqBAIJ,0, 499b6490206SBarry Smith 0,0, 500de6a44a3SBarry Smith MatLUFactor_SeqBAIJ,0, 501b6490206SBarry Smith 0, 502b6490206SBarry Smith 0, 503b6490206SBarry Smith MatGetInfo_SeqBAIJ,0, 504b6490206SBarry Smith 0,0,MatNorm_SeqBAIJ, 505b6490206SBarry Smith 0,0, 506b6490206SBarry Smith 0, 507b6490206SBarry Smith MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0, 508b6490206SBarry Smith MatGetReordering_SeqBAIJ, 5097fc0212eSBarry Smith MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 510b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 511de6a44a3SBarry Smith MatILUFactorSymbolic_SeqBAIJ,0, 512b6490206SBarry Smith 0,0,/*MatConvert_SeqBAIJ*/ 0, 513b6490206SBarry Smith 0,0, 514b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 515b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 516b6490206SBarry Smith 0,0, 517b6490206SBarry Smith 0,0, 518b6490206SBarry Smith MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 519b6490206SBarry Smith 0}; 5202593348eSBarry Smith 5212593348eSBarry Smith /*@C 522b6490206SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format 5232593348eSBarry Smith (the default parallel PETSc format). For good matrix assembly performance 5242593348eSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 5252593348eSBarry Smith (or nzz). By setting these parameters accurately, performance can be 5262593348eSBarry Smith increased by more than a factor of 50. 5272593348eSBarry Smith 5282593348eSBarry Smith Input Parameters: 5292593348eSBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 530b6490206SBarry Smith . bs - size of block 5312593348eSBarry Smith . m - number of rows 5322593348eSBarry Smith . n - number of columns 533b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 534b6490206SBarry Smith . nzz - number of block nonzeros per block row or PETSC_NULL 535b6490206SBarry Smith (possibly different for each row) 5362593348eSBarry Smith 5372593348eSBarry Smith Output Parameter: 5382593348eSBarry Smith . A - the matrix 5392593348eSBarry Smith 5402593348eSBarry Smith Notes: 541b6490206SBarry Smith The BAIJ format, is fully compatible with standard Fortran 77 5422593348eSBarry Smith storage. That is, the stored row and column indices can begin at 5432593348eSBarry Smith either one (as in Fortran) or zero. See the users manual for details. 5442593348eSBarry Smith 5452593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 5462593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 5472593348eSBarry Smith allocation. For additional details, see the users manual chapter on 5482593348eSBarry Smith matrices and the file $(PETSC_DIR)/Performance. 5492593348eSBarry Smith 5502593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 5512593348eSBarry Smith @*/ 552b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 5532593348eSBarry Smith { 5542593348eSBarry Smith Mat B; 555b6490206SBarry Smith Mat_SeqBAIJ *b; 556b6490206SBarry Smith int i,len,ierr, flg,mbs = m/bs; 557b6490206SBarry Smith 558b6490206SBarry Smith if (mbs*bs != m) 559b6490206SBarry Smith SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize"); 5602593348eSBarry Smith 5612593348eSBarry Smith *A = 0; 562b6490206SBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm); 5632593348eSBarry Smith PLogObjectCreate(B); 564b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 5652593348eSBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 5667fc0212eSBarry Smith switch (bs) { 5677fc0212eSBarry Smith case 1: 5687fc0212eSBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 5697fc0212eSBarry Smith break; 5704eeb42bcSBarry Smith case 2: 5714eeb42bcSBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 5726d84be18SBarry Smith break; 573f327dd97SBarry Smith case 3: 574f327dd97SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 5754eeb42bcSBarry Smith break; 5766d84be18SBarry Smith case 4: 5776d84be18SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 5786d84be18SBarry Smith break; 5796d84be18SBarry Smith case 5: 5806d84be18SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 5816d84be18SBarry Smith break; 5827fc0212eSBarry Smith } 583b6490206SBarry Smith B->destroy = MatDestroy_SeqBAIJ; 584b6490206SBarry Smith B->view = MatView_SeqBAIJ; 5852593348eSBarry Smith B->factor = 0; 5862593348eSBarry Smith B->lupivotthreshold = 1.0; 5872593348eSBarry Smith b->row = 0; 5882593348eSBarry Smith b->col = 0; 5892593348eSBarry Smith b->reallocs = 0; 5902593348eSBarry Smith 5912593348eSBarry Smith b->m = m; 592b6490206SBarry Smith b->mbs = mbs; 5932593348eSBarry Smith b->n = n; 594b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 5952593348eSBarry Smith if (nnz == PETSC_NULL) { 596de6a44a3SBarry Smith if (nz == PETSC_DEFAULT) nz = 5; 5972593348eSBarry Smith else if (nz <= 0) nz = 1; 598b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 599b6490206SBarry Smith nz = nz*mbs; 6002593348eSBarry Smith } 6012593348eSBarry Smith else { 6022593348eSBarry Smith nz = 0; 603b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 6042593348eSBarry Smith } 6052593348eSBarry Smith 6062593348eSBarry Smith /* allocate the matrix space */ 607b6490206SBarry Smith len = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int); 6082593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 609b6490206SBarry Smith PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar)); 610b6490206SBarry Smith b->j = (int *) (b->a + nz*bs*bs); 6112593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 6122593348eSBarry Smith b->i = b->j + nz; 6132593348eSBarry Smith b->singlemalloc = 1; 6142593348eSBarry Smith 615de6a44a3SBarry Smith b->i[0] = 0; 616b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 6172593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 6182593348eSBarry Smith } 6192593348eSBarry Smith 620b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 621b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 622b6490206SBarry Smith PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 623b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 6242593348eSBarry Smith 625b6490206SBarry Smith b->bs = bs; 626b6490206SBarry Smith b->mbs = mbs; 6272593348eSBarry Smith b->nz = 0; 6282593348eSBarry Smith b->maxnz = nz; 6292593348eSBarry Smith b->sorted = 0; 6302593348eSBarry Smith b->roworiented = 1; 6312593348eSBarry Smith b->nonew = 0; 6322593348eSBarry Smith b->diag = 0; 6332593348eSBarry Smith b->solve_work = 0; 634de6a44a3SBarry Smith b->mult_work = 0; 6352593348eSBarry Smith b->spptr = 0; 6362593348eSBarry Smith 6372593348eSBarry Smith *A = B; 6382593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 6392593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 6402593348eSBarry Smith return 0; 6412593348eSBarry Smith } 6422593348eSBarry Smith 643b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 6442593348eSBarry Smith { 6452593348eSBarry Smith Mat C; 646b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 647de6a44a3SBarry Smith int i,len, mbs = a->mbs, bs = a->bs,nz = a->nz; 648de6a44a3SBarry Smith 649de6a44a3SBarry Smith if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix"); 6502593348eSBarry Smith 6512593348eSBarry Smith *B = 0; 652b6490206SBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm); 6532593348eSBarry Smith PLogObjectCreate(C); 654b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 6552593348eSBarry Smith PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 656b6490206SBarry Smith C->destroy = MatDestroy_SeqBAIJ; 657b6490206SBarry Smith C->view = MatView_SeqBAIJ; 6582593348eSBarry Smith C->factor = A->factor; 6592593348eSBarry Smith c->row = 0; 6602593348eSBarry Smith c->col = 0; 6612593348eSBarry Smith C->assembled = PETSC_TRUE; 6622593348eSBarry Smith 6632593348eSBarry Smith c->m = a->m; 6642593348eSBarry Smith c->n = a->n; 665b6490206SBarry Smith c->bs = a->bs; 666b6490206SBarry Smith c->mbs = a->mbs; 667b6490206SBarry Smith c->nbs = a->nbs; 6682593348eSBarry Smith 669b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 670b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 671b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 6722593348eSBarry Smith c->imax[i] = a->imax[i]; 6732593348eSBarry Smith c->ilen[i] = a->ilen[i]; 6742593348eSBarry Smith } 6752593348eSBarry Smith 6762593348eSBarry Smith /* allocate the matrix space */ 6772593348eSBarry Smith c->singlemalloc = 1; 678de6a44a3SBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int)); 6792593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 680de6a44a3SBarry Smith c->j = (int *) (c->a + nz*bs*bs); 681de6a44a3SBarry Smith c->i = c->j + nz; 682b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 683b6490206SBarry Smith if (mbs > 0) { 684de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 6852593348eSBarry Smith if (cpvalues == COPY_VALUES) { 686de6a44a3SBarry Smith PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar)); 6872593348eSBarry Smith } 6882593348eSBarry Smith } 6892593348eSBarry Smith 690b6490206SBarry Smith PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 6912593348eSBarry Smith c->sorted = a->sorted; 6922593348eSBarry Smith c->roworiented = a->roworiented; 6932593348eSBarry Smith c->nonew = a->nonew; 6942593348eSBarry Smith 6952593348eSBarry Smith if (a->diag) { 696b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 697b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 698b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 6992593348eSBarry Smith c->diag[i] = a->diag[i]; 7002593348eSBarry Smith } 7012593348eSBarry Smith } 7022593348eSBarry Smith else c->diag = 0; 7032593348eSBarry Smith c->nz = a->nz; 7042593348eSBarry Smith c->maxnz = a->maxnz; 7052593348eSBarry Smith c->solve_work = 0; 7062593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 7077fc0212eSBarry Smith c->mult_work = 0; 7082593348eSBarry Smith *B = C; 7092593348eSBarry Smith return 0; 7102593348eSBarry Smith } 7112593348eSBarry Smith 712b6490206SBarry Smith int MatLoad_SeqBAIJ(Viewer bview,MatType type,Mat *A) 7132593348eSBarry Smith { 714b6490206SBarry Smith Mat_SeqBAIJ *a; 7152593348eSBarry Smith Mat B; 716de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 717b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 71835aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 719de6a44a3SBarry Smith int *masked, nmask,tmp,ishift,bs2; 720b6490206SBarry Smith Scalar *aa; 7212593348eSBarry Smith PetscObject vobj = (PetscObject) bview; 7222593348eSBarry Smith MPI_Comm comm = vobj->comm; 7232593348eSBarry Smith 72435aab85fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 725de6a44a3SBarry Smith bs2 = bs*bs; 726b6490206SBarry Smith 7272593348eSBarry Smith MPI_Comm_size(comm,&size); 728b6490206SBarry Smith if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor"); 7292593348eSBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 7302593348eSBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 731de6a44a3SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object"); 7322593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 7332593348eSBarry Smith 734b6490206SBarry Smith if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices"); 73535aab85fSBarry Smith 73635aab85fSBarry Smith /* 73735aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 73835aab85fSBarry Smith divisible by the blocksize 73935aab85fSBarry Smith */ 740b6490206SBarry Smith mbs = M/bs; 74135aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 74235aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 74335aab85fSBarry Smith else mbs++; 74435aab85fSBarry Smith if (extra_rows) { 74535aab85fSBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize"); 74635aab85fSBarry Smith } 747b6490206SBarry Smith 7482593348eSBarry Smith /* read in row lengths */ 74935aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 7502593348eSBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 75135aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 7522593348eSBarry Smith 753b6490206SBarry Smith /* read in column indices */ 75435aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 755b6490206SBarry Smith ierr = SYRead(fd,jj,nz,SYINT); CHKERRQ(ierr); 75635aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 757b6490206SBarry Smith 758b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 759b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 760b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 76135aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 76235aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 76335aab85fSBarry Smith masked = mask + mbs; 764b6490206SBarry Smith rowcount = 0; nzcount = 0; 765b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 76635aab85fSBarry Smith nmask = 0; 767b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 768b6490206SBarry Smith kmax = rowlengths[rowcount]; 769b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 77035aab85fSBarry Smith tmp = jj[nzcount++]/bs; 77135aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 772b6490206SBarry Smith } 773b6490206SBarry Smith rowcount++; 774b6490206SBarry Smith } 77535aab85fSBarry Smith browlengths[i] += nmask; 77635aab85fSBarry Smith /* zero out the mask elements we set */ 77735aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 778b6490206SBarry Smith } 779b6490206SBarry Smith 7802593348eSBarry Smith /* create our matrix */ 78135aab85fSBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A); 78235aab85fSBarry Smith CHKERRQ(ierr); 7832593348eSBarry Smith B = *A; 784b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 7852593348eSBarry Smith 7862593348eSBarry Smith /* set matrix "i" values */ 787de6a44a3SBarry Smith a->i[0] = 0; 788b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 789b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 790b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 7912593348eSBarry Smith } 792b6490206SBarry Smith a->nz = 0; 793b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 7942593348eSBarry Smith 795b6490206SBarry Smith /* read in nonzero values */ 79635aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 797b6490206SBarry Smith ierr = SYRead(fd,aa,nz,SYSCALAR); CHKERRQ(ierr); 79835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 799b6490206SBarry Smith 800b6490206SBarry Smith /* set "a" and "j" values into matrix */ 801b6490206SBarry Smith nzcount = 0; jcount = 0; 802b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 803b6490206SBarry Smith nzcountb = nzcount; 80435aab85fSBarry Smith nmask = 0; 805b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 806b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 807b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 80835aab85fSBarry Smith tmp = jj[nzcount++]/bs; 80935aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 810b6490206SBarry Smith } 811b6490206SBarry Smith rowcount++; 812b6490206SBarry Smith } 813de6a44a3SBarry Smith /* sort the masked values */ 814de6a44a3SBarry Smith SYIsort(nmask,masked); 815de6a44a3SBarry Smith 816b6490206SBarry Smith /* set "j" values into matrix */ 817b6490206SBarry Smith maskcount = 1; 81835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 81935aab85fSBarry Smith a->j[jcount++] = masked[j]; 820de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 821b6490206SBarry Smith } 822b6490206SBarry Smith /* set "a" values into matrix */ 823de6a44a3SBarry Smith ishift = bs2*a->i[i]; 824b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 825b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 826b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 827de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 828de6a44a3SBarry Smith block = mask[tmp] - 1; 829de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 830de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 831b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 832b6490206SBarry Smith } 833b6490206SBarry Smith } 83435aab85fSBarry Smith /* zero out the mask elements we set */ 83535aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 836b6490206SBarry Smith } 83735aab85fSBarry Smith if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix"); 838b6490206SBarry Smith 839b6490206SBarry Smith PetscFree(rowlengths); 840b6490206SBarry Smith PetscFree(browlengths); 841b6490206SBarry Smith PetscFree(aa); 842b6490206SBarry Smith PetscFree(jj); 843b6490206SBarry Smith PetscFree(mask); 844b6490206SBarry Smith 845b6490206SBarry Smith B->assembled = PETSC_TRUE; 846de6a44a3SBarry Smith 847de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 848de6a44a3SBarry Smith if (flg) { 849de6a44a3SBarry Smith Viewer viewer; 850de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 851de6a44a3SBarry Smith ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO,0);CHKERRQ(ierr); 852de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 853de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 854de6a44a3SBarry Smith } 855de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 856de6a44a3SBarry Smith if (flg) { 857de6a44a3SBarry Smith Viewer viewer; 858de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 859de6a44a3SBarry Smith ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr); 860de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 861de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 862de6a44a3SBarry Smith } 863de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr); 864de6a44a3SBarry Smith if (flg) { 865de6a44a3SBarry Smith Viewer viewer; 866de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 867de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 868de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 869de6a44a3SBarry Smith } 870de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr); 871de6a44a3SBarry Smith if (flg) { 872de6a44a3SBarry Smith Viewer viewer; 873de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 874de6a44a3SBarry Smith ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_MATLAB,"M");CHKERRQ(ierr); 875de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 876de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 877de6a44a3SBarry Smith } 878de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 879de6a44a3SBarry Smith if (flg) { 880de6a44a3SBarry Smith Draw win; 881de6a44a3SBarry Smith ierr = DrawOpenX(B->comm,0,0,0,0,300,300,&win); CHKERRQ(ierr); 882de6a44a3SBarry Smith ierr = MatView(B,(Viewer)win); CHKERRQ(ierr); 883de6a44a3SBarry Smith ierr = DrawSyncFlush(win); CHKERRQ(ierr); 884de6a44a3SBarry Smith ierr = DrawDestroy(win); CHKERRQ(ierr); 885de6a44a3SBarry Smith } 8862593348eSBarry Smith return 0; 8872593348eSBarry Smith } 8882593348eSBarry Smith 8892593348eSBarry Smith 8902593348eSBarry Smith 891