1b6490206SBarry Smith 22593348eSBarry Smith #ifndef lint 3*de6a44a3SBarry Smith static char vcid[] = "$Id: baij.c,v 1.3 1996/02/13 23:29:47 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 15b6490206SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(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*de6a44a3SBarry Smith int ierr, *ia, *ja,n = a->mbs,*idx,i; 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*de6a44a3SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,&ia,&ja);CHKERRQ(ierr); 41*de6a44a3SBarry Smith ierr = MatGetReordering_IJ(n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 422593348eSBarry Smith 432593348eSBarry Smith PetscFree(ia); PetscFree(ja); 442593348eSBarry Smith return 0; 452593348eSBarry Smith } 462593348eSBarry Smith 47*de6a44a3SBarry Smith /* 48*de6a44a3SBarry Smith Adds diagonal pointers to sparse matrix structure. 49*de6a44a3SBarry Smith */ 50*de6a44a3SBarry Smith 51*de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 52*de6a44a3SBarry Smith { 53*de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 54*de6a44a3SBarry Smith int i,j, *diag, m = a->m; 55*de6a44a3SBarry Smith 56*de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 57*de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 58*de6a44a3SBarry Smith for ( i=0; i<a->m; i++ ) { 59*de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 60*de6a44a3SBarry Smith if (a->j[j] == i) { 61*de6a44a3SBarry Smith diag[i] = j; 62*de6a44a3SBarry Smith break; 63*de6a44a3SBarry Smith } 64*de6a44a3SBarry Smith } 65*de6a44a3SBarry Smith } 66*de6a44a3SBarry Smith a->diag = diag; 67*de6a44a3SBarry Smith return 0; 68*de6a44a3SBarry Smith } 692593348eSBarry Smith 702593348eSBarry Smith #include "draw.h" 712593348eSBarry Smith #include "pinclude/pviewer.h" 722593348eSBarry Smith #include "sysio.h" 732593348eSBarry Smith 74b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 752593348eSBarry Smith { 76b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 77b6490206SBarry Smith int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l; 78b6490206SBarry Smith Scalar *aa; 792593348eSBarry Smith 802593348eSBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 812593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 822593348eSBarry Smith col_lens[0] = MAT_COOKIE; 832593348eSBarry Smith col_lens[1] = a->m; 842593348eSBarry Smith col_lens[2] = a->n; 85b6490206SBarry Smith col_lens[3] = a->nz*bs*bs; 862593348eSBarry Smith 872593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 88b6490206SBarry Smith count = 0; 89b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 90b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 91b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 92b6490206SBarry Smith } 932593348eSBarry Smith } 942593348eSBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 952593348eSBarry Smith PetscFree(col_lens); 962593348eSBarry Smith 972593348eSBarry Smith /* store column indices (zero start index) */ 98b6490206SBarry Smith jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj); 99b6490206SBarry Smith count = 0; 100b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 101b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 102b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 103b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 104b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 1052593348eSBarry Smith } 1062593348eSBarry Smith } 107b6490206SBarry Smith } 108b6490206SBarry Smith } 109b6490206SBarry Smith ierr = SYWrite(fd,jj,bs*bs*a->nz,SYINT,0); CHKERRQ(ierr); 110b6490206SBarry Smith PetscFree(jj); 1112593348eSBarry Smith 1122593348eSBarry Smith /* store nonzero values */ 113b6490206SBarry Smith aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa); 114b6490206SBarry Smith count = 0; 115b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 116b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 117b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 118b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 119b6490206SBarry Smith aa[count++] = a->a[bs*bs*k + l*bs + j]; 120b6490206SBarry Smith } 121b6490206SBarry Smith } 122b6490206SBarry Smith } 123b6490206SBarry Smith } 124b6490206SBarry Smith ierr = SYWrite(fd,aa,bs*bs*a->nz,SYSCALAR,0); CHKERRQ(ierr); 125b6490206SBarry Smith PetscFree(aa); 1262593348eSBarry Smith return 0; 1272593348eSBarry Smith } 1282593348eSBarry Smith 129b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 1302593348eSBarry Smith { 131b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 132*de6a44a3SBarry Smith int ierr, i,j,format,bs = a->bs,k,l; 1332593348eSBarry Smith FILE *fd; 1342593348eSBarry Smith char *outputname; 1352593348eSBarry Smith 1362593348eSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 1372593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 1382593348eSBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 1392593348eSBarry Smith if (format == FILE_FORMAT_INFO) { 1402593348eSBarry Smith /* no need to print additional information */ ; 1412593348eSBarry Smith } 1422593348eSBarry Smith else if (format == FILE_FORMAT_MATLAB) { 143b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported"); 1442593348eSBarry Smith } 1452593348eSBarry Smith else { 146b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 147b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 148b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 149b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 150b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 151*de6a44a3SBarry Smith fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs*bs*k + l*bs + j]); 1522593348eSBarry Smith } 1532593348eSBarry Smith } 1542593348eSBarry Smith fprintf(fd,"\n"); 1552593348eSBarry Smith } 1562593348eSBarry Smith } 157b6490206SBarry Smith } 1582593348eSBarry Smith fflush(fd); 1592593348eSBarry Smith return 0; 1602593348eSBarry Smith } 1612593348eSBarry Smith 162b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer) 1632593348eSBarry Smith { 1642593348eSBarry Smith Mat A = (Mat) obj; 1652593348eSBarry Smith PetscObject vobj = (PetscObject) viewer; 1662593348eSBarry Smith 1672593348eSBarry Smith if (!viewer) { 1682593348eSBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 1692593348eSBarry Smith } 1702593348eSBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 1712593348eSBarry Smith if (vobj->type == MATLAB_VIEWER) { 172b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported"); 1732593348eSBarry Smith } 1742593348eSBarry Smith else if (vobj->type == ASCII_FILE_VIEWER||vobj->type == ASCII_FILES_VIEWER){ 175b6490206SBarry Smith return MatView_SeqBAIJ_ASCII(A,viewer); 1762593348eSBarry Smith } 1772593348eSBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 178b6490206SBarry Smith return MatView_SeqBAIJ_Binary(A,viewer); 1792593348eSBarry Smith } 1802593348eSBarry Smith } 1812593348eSBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 1822593348eSBarry Smith if (vobj->type == NULLWINDOW) return 0; 183b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported"); 1842593348eSBarry Smith } 1852593348eSBarry Smith return 0; 1862593348eSBarry Smith } 187b6490206SBarry Smith 188b6490206SBarry Smith 189b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A) 1902593348eSBarry Smith { 191b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 192*de6a44a3SBarry Smith PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar)); 1932593348eSBarry Smith return 0; 1942593348eSBarry Smith } 1952593348eSBarry Smith 196b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj) 1972593348eSBarry Smith { 1982593348eSBarry Smith Mat A = (Mat) obj; 199b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2002593348eSBarry Smith 2012593348eSBarry Smith #if defined(PETSC_LOG) 2022593348eSBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 2032593348eSBarry Smith #endif 2042593348eSBarry Smith PetscFree(a->a); 2052593348eSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 2062593348eSBarry Smith if (a->diag) PetscFree(a->diag); 2072593348eSBarry Smith if (a->ilen) PetscFree(a->ilen); 2082593348eSBarry Smith if (a->imax) PetscFree(a->imax); 2092593348eSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 210*de6a44a3SBarry Smith if (a->mult_work) PetscFree(a->mult_work); 2112593348eSBarry Smith PetscFree(a); 2122593348eSBarry Smith PLogObjectDestroy(A); 2132593348eSBarry Smith PetscHeaderDestroy(A); 2142593348eSBarry Smith return 0; 2152593348eSBarry Smith } 2162593348eSBarry Smith 217b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op) 2182593348eSBarry Smith { 219b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2202593348eSBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 2212593348eSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 2222593348eSBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 2232593348eSBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 2242593348eSBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 2252593348eSBarry Smith else if (op == ROWS_SORTED || 2262593348eSBarry Smith op == SYMMETRIC_MATRIX || 2272593348eSBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 2282593348eSBarry Smith op == YES_NEW_DIAGONALS) 229b6490206SBarry Smith PLogInfo((PetscObject)A,"Info:MatSetOption_SeqBAIJ:Option ignored\n"); 2302593348eSBarry Smith else if (op == NO_NEW_DIAGONALS) 231b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");} 2322593348eSBarry Smith else 233b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");} 2342593348eSBarry Smith return 0; 2352593348eSBarry Smith } 2362593348eSBarry Smith 2372593348eSBarry Smith 2382593348eSBarry Smith /* -------------------------------------------------------*/ 2392593348eSBarry Smith /* Should check that shapes of vectors and matrices match */ 2402593348eSBarry Smith /* -------------------------------------------------------*/ 241b6490206SBarry Smith #include "pinclude/plapack.h" 242b6490206SBarry Smith 243b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy) 2442593348eSBarry Smith { 245b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 246b6490206SBarry Smith Scalar *xg,*yg; 247b6490206SBarry Smith register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5; 248b6490206SBarry Smith register Scalar x1,x2,x3,x4,x5; 249b6490206SBarry Smith int mbs = a->mbs, m = a->m, i, *idx,*ii; 250b6490206SBarry Smith int bs = a->bs,j,n,bs2 = bs*bs; 2512593348eSBarry Smith 252b6490206SBarry Smith VecGetArray(xx,&xg); x = xg; VecGetArray(yy,&yg); y = yg; 253b6490206SBarry Smith PetscMemzero(y,m*sizeof(Scalar)); 254b6490206SBarry Smith x = x; 2552593348eSBarry Smith idx = a->j; 2562593348eSBarry Smith v = a->a; 2572593348eSBarry Smith ii = a->i; 258b6490206SBarry Smith 259b6490206SBarry Smith switch (bs) { 260b6490206SBarry Smith case 1: 2612593348eSBarry Smith for ( i=0; i<m; i++ ) { 2622593348eSBarry Smith n = ii[1] - ii[0]; ii++; 2632593348eSBarry Smith sum = 0.0; 2642593348eSBarry Smith while (n--) sum += *v++ * x[*idx++]; 2652593348eSBarry Smith y[i] = sum; 2662593348eSBarry Smith } 2672593348eSBarry Smith break; 268b6490206SBarry Smith case 2: 269b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 270*de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 271b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; 272b6490206SBarry Smith for ( j=0; j<n; j++ ) { 273b6490206SBarry Smith xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 274b6490206SBarry Smith sum1 += v[0]*x1 + v[2]*x2; 275b6490206SBarry Smith sum2 += v[1]*x1 + v[3]*x2; 276b6490206SBarry Smith v += 4; 277b6490206SBarry Smith } 278b6490206SBarry Smith y[0] += sum1; y[1] += sum2; 279b6490206SBarry Smith y += 2; 280b6490206SBarry Smith } 281b6490206SBarry Smith break; 282b6490206SBarry Smith case 3: 283b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 284*de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 285b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 286b6490206SBarry Smith for ( j=0; j<n; j++ ) { 287b6490206SBarry Smith xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 288b6490206SBarry Smith sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 289b6490206SBarry Smith sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 290b6490206SBarry Smith sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 291b6490206SBarry Smith v += 9; 292b6490206SBarry Smith } 293b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; 294b6490206SBarry Smith y += 3; 295b6490206SBarry Smith } 296b6490206SBarry Smith break; 297b6490206SBarry Smith case 4: 298b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 299*de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 300b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 301b6490206SBarry Smith for ( j=0; j<n; j++ ) { 302b6490206SBarry Smith xb = x + 4*(*idx++); 303b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 304b6490206SBarry Smith sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 305b6490206SBarry Smith sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 306b6490206SBarry Smith sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 307b6490206SBarry Smith sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 308b6490206SBarry Smith v += 16; 309b6490206SBarry Smith } 310b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; 311b6490206SBarry Smith y += 4; 312b6490206SBarry Smith } 313b6490206SBarry Smith break; 314b6490206SBarry Smith case 5: 315b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 316*de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 317b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 318b6490206SBarry Smith for ( j=0; j<n; j++ ) { 319b6490206SBarry Smith xb = x + 5*(*idx++); 320b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 321b6490206SBarry Smith sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 322b6490206SBarry Smith sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 323b6490206SBarry Smith sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 324b6490206SBarry Smith sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 325b6490206SBarry Smith sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 326b6490206SBarry Smith v += 25; 327b6490206SBarry Smith } 328b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5; 329b6490206SBarry Smith y += 5; 330b6490206SBarry Smith } 331b6490206SBarry Smith break; 332b6490206SBarry Smith /* block sizes larger then 5 by 5 are handled by BLAS */ 333b6490206SBarry Smith default: { 334*de6a44a3SBarry Smith int _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt; 335*de6a44a3SBarry Smith if (!a->mult_work) { 336*de6a44a3SBarry Smith a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar)); 337*de6a44a3SBarry Smith CHKPTRQ(a->mult_work); 338*de6a44a3SBarry Smith } 339*de6a44a3SBarry Smith work = a->mult_work; 340b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 341*de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 342*de6a44a3SBarry Smith ncols = n*bs; 343*de6a44a3SBarry Smith workt = work; 344b6490206SBarry Smith for ( j=0; j<n; j++ ) { 345b6490206SBarry Smith xb = x + bs*(*idx++); 346*de6a44a3SBarry Smith for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 347*de6a44a3SBarry Smith workt += bs; 348b6490206SBarry Smith } 349*de6a44a3SBarry Smith LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One); 350*de6a44a3SBarry Smith v += n*bs2; 351b6490206SBarry Smith y += bs; 3522593348eSBarry Smith } 3532593348eSBarry Smith } 3542593348eSBarry Smith } 355b6490206SBarry Smith PLogFlops(2*bs2*a->nz - m); 3562593348eSBarry Smith return 0; 3572593348eSBarry Smith } 3582593348eSBarry Smith 359*de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem) 3602593348eSBarry Smith { 361b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 36235aab85fSBarry Smith *nz = a->bs*a->bs*a->nz; 363*de6a44a3SBarry Smith *nza = a->maxnz; 3642593348eSBarry Smith *mem = (int)A->mem; 3652593348eSBarry Smith return 0; 3662593348eSBarry Smith } 3672593348eSBarry Smith 368b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 369b6490206SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ(Mat,Mat*); 370b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 371b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec); 3722593348eSBarry Smith 373b6490206SBarry Smith static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 3742593348eSBarry Smith { 375b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3762593348eSBarry Smith *m = a->m; *n = a->n; 3772593348eSBarry Smith return 0; 3782593348eSBarry Smith } 3792593348eSBarry Smith 380b6490206SBarry Smith static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 3812593348eSBarry Smith { 382b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3832593348eSBarry Smith *m = 0; *n = a->m; 3842593348eSBarry Smith return 0; 3852593348eSBarry Smith } 386b6490206SBarry Smith 387b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 3882593348eSBarry Smith { 389b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3902593348eSBarry Smith Scalar *v = a->a; 3912593348eSBarry Smith double sum = 0.0; 392b6490206SBarry Smith int i; 3932593348eSBarry Smith 3942593348eSBarry Smith if (type == NORM_FROBENIUS) { 3952593348eSBarry Smith for (i=0; i<a->nz; i++ ) { 3962593348eSBarry Smith #if defined(PETSC_COMPLEX) 3972593348eSBarry Smith sum += real(conj(*v)*(*v)); v++; 3982593348eSBarry Smith #else 3992593348eSBarry Smith sum += (*v)*(*v); v++; 4002593348eSBarry Smith #endif 4012593348eSBarry Smith } 4022593348eSBarry Smith *norm = sqrt(sum); 4032593348eSBarry Smith } 4042593348eSBarry Smith else { 405b6490206SBarry Smith SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet"); 4062593348eSBarry Smith } 4072593348eSBarry Smith return 0; 4082593348eSBarry Smith } 4092593348eSBarry Smith 4102593348eSBarry Smith /* 4112593348eSBarry Smith note: This can only work for identity for row and col. It would 4122593348eSBarry Smith be good to check this and otherwise generate an error. 4132593348eSBarry Smith */ 414b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 4152593348eSBarry Smith { 416b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 4172593348eSBarry Smith Mat outA; 418*de6a44a3SBarry Smith int ierr; 4192593348eSBarry Smith 420b6490206SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported"); 4212593348eSBarry Smith 4222593348eSBarry Smith outA = inA; 4232593348eSBarry Smith inA->factor = FACTOR_LU; 4242593348eSBarry Smith a->row = row; 4252593348eSBarry Smith a->col = col; 4262593348eSBarry Smith 427*de6a44a3SBarry Smith a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 4282593348eSBarry Smith 4292593348eSBarry Smith if (!a->diag) { 430b6490206SBarry Smith ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 4312593348eSBarry Smith } 432b6490206SBarry Smith ierr = MatLUFactorNumeric_SeqBAIJ(inA,&outA); CHKERRQ(ierr); 4332593348eSBarry Smith return 0; 4342593348eSBarry Smith } 4352593348eSBarry Smith 436b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 4372593348eSBarry Smith { 438b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 439b6490206SBarry Smith int one = 1, totalnz = a->bs*a->bs*a->nz; 440b6490206SBarry Smith BLscal_( &totalnz, alpha, a->a, &one ); 441b6490206SBarry Smith PLogFlops(totalnz); 4422593348eSBarry Smith return 0; 4432593348eSBarry Smith } 4442593348eSBarry Smith 445b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A) 4462593348eSBarry Smith { 4472593348eSBarry Smith static int called = 0; 4482593348eSBarry Smith MPI_Comm comm = A->comm; 4492593348eSBarry Smith 4502593348eSBarry Smith if (called) return 0; else called = 1; 4512593348eSBarry Smith return 0; 4522593348eSBarry Smith } 4532593348eSBarry Smith /* -------------------------------------------------------------------*/ 454b6490206SBarry Smith static struct _MatOps MatOps = {0, 455b6490206SBarry Smith 0,0, 456b6490206SBarry Smith MatMult_SeqBAIJ,0, 457b6490206SBarry Smith 0,0, 458*de6a44a3SBarry Smith MatSolve_SeqBAIJ,0, 459b6490206SBarry Smith 0,0, 460*de6a44a3SBarry Smith MatLUFactor_SeqBAIJ,0, 461b6490206SBarry Smith 0, 462b6490206SBarry Smith 0, 463b6490206SBarry Smith MatGetInfo_SeqBAIJ,0, 464b6490206SBarry Smith 0,0,MatNorm_SeqBAIJ, 465b6490206SBarry Smith 0,0, 466b6490206SBarry Smith 0, 467b6490206SBarry Smith MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0, 468b6490206SBarry Smith MatGetReordering_SeqBAIJ, 469*de6a44a3SBarry Smith MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ,0,0, 470b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 471*de6a44a3SBarry Smith MatILUFactorSymbolic_SeqBAIJ,0, 472b6490206SBarry Smith 0,0,/*MatConvert_SeqBAIJ*/ 0, 473b6490206SBarry Smith 0,0, 474b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 475b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 476b6490206SBarry Smith 0,0, 477b6490206SBarry Smith 0,0, 478b6490206SBarry Smith MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 479b6490206SBarry Smith 0}; 4802593348eSBarry Smith 4812593348eSBarry Smith /*@C 482b6490206SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format 4832593348eSBarry Smith (the default parallel PETSc format). For good matrix assembly performance 4842593348eSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 4852593348eSBarry Smith (or nzz). By setting these parameters accurately, performance can be 4862593348eSBarry Smith increased by more than a factor of 50. 4872593348eSBarry Smith 4882593348eSBarry Smith Input Parameters: 4892593348eSBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 490b6490206SBarry Smith . bs - size of block 4912593348eSBarry Smith . m - number of rows 4922593348eSBarry Smith . n - number of columns 493b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 494b6490206SBarry Smith . nzz - number of block nonzeros per block row or PETSC_NULL 495b6490206SBarry Smith (possibly different for each row) 4962593348eSBarry Smith 4972593348eSBarry Smith Output Parameter: 4982593348eSBarry Smith . A - the matrix 4992593348eSBarry Smith 5002593348eSBarry Smith Notes: 501b6490206SBarry Smith The BAIJ format, is fully compatible with standard Fortran 77 5022593348eSBarry Smith storage. That is, the stored row and column indices can begin at 5032593348eSBarry Smith either one (as in Fortran) or zero. See the users manual for details. 5042593348eSBarry Smith 5052593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 5062593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 5072593348eSBarry Smith allocation. For additional details, see the users manual chapter on 5082593348eSBarry Smith matrices and the file $(PETSC_DIR)/Performance. 5092593348eSBarry Smith 5102593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 5112593348eSBarry Smith @*/ 512b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 5132593348eSBarry Smith { 5142593348eSBarry Smith Mat B; 515b6490206SBarry Smith Mat_SeqBAIJ *b; 516b6490206SBarry Smith int i,len,ierr, flg,mbs = m/bs; 517b6490206SBarry Smith 518b6490206SBarry Smith if (mbs*bs != m) 519b6490206SBarry Smith SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize"); 5202593348eSBarry Smith 5212593348eSBarry Smith *A = 0; 522b6490206SBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm); 5232593348eSBarry Smith PLogObjectCreate(B); 524b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 5252593348eSBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 526b6490206SBarry Smith B->destroy = MatDestroy_SeqBAIJ; 527b6490206SBarry Smith B->view = MatView_SeqBAIJ; 5282593348eSBarry Smith B->factor = 0; 5292593348eSBarry Smith B->lupivotthreshold = 1.0; 5302593348eSBarry Smith b->row = 0; 5312593348eSBarry Smith b->col = 0; 5322593348eSBarry Smith b->reallocs = 0; 5332593348eSBarry Smith 5342593348eSBarry Smith b->m = m; 535b6490206SBarry Smith b->mbs = mbs; 5362593348eSBarry Smith b->n = n; 537b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 5382593348eSBarry Smith if (nnz == PETSC_NULL) { 539*de6a44a3SBarry Smith if (nz == PETSC_DEFAULT) nz = 5; 5402593348eSBarry Smith else if (nz <= 0) nz = 1; 541b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 542b6490206SBarry Smith nz = nz*mbs; 5432593348eSBarry Smith } 5442593348eSBarry Smith else { 5452593348eSBarry Smith nz = 0; 546b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 5472593348eSBarry Smith } 5482593348eSBarry Smith 5492593348eSBarry Smith /* allocate the matrix space */ 550b6490206SBarry Smith len = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int); 5512593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 552b6490206SBarry Smith PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar)); 553b6490206SBarry Smith b->j = (int *) (b->a + nz*bs*bs); 5542593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 5552593348eSBarry Smith b->i = b->j + nz; 5562593348eSBarry Smith b->singlemalloc = 1; 5572593348eSBarry Smith 558*de6a44a3SBarry Smith b->i[0] = 0; 559b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 5602593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 5612593348eSBarry Smith } 5622593348eSBarry Smith 563b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 564b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 565b6490206SBarry Smith PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 566b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 5672593348eSBarry Smith 568b6490206SBarry Smith b->bs = bs; 569b6490206SBarry Smith b->mbs = mbs; 5702593348eSBarry Smith b->nz = 0; 5712593348eSBarry Smith b->maxnz = nz; 5722593348eSBarry Smith b->sorted = 0; 5732593348eSBarry Smith b->roworiented = 1; 5742593348eSBarry Smith b->nonew = 0; 5752593348eSBarry Smith b->diag = 0; 5762593348eSBarry Smith b->solve_work = 0; 577*de6a44a3SBarry Smith b->mult_work = 0; 5782593348eSBarry Smith b->spptr = 0; 5792593348eSBarry Smith 5802593348eSBarry Smith *A = B; 5812593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 5822593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 5832593348eSBarry Smith return 0; 5842593348eSBarry Smith } 5852593348eSBarry Smith 586b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 5872593348eSBarry Smith { 5882593348eSBarry Smith Mat C; 589b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 590*de6a44a3SBarry Smith int i,len, mbs = a->mbs, bs = a->bs,nz = a->nz; 591*de6a44a3SBarry Smith 592*de6a44a3SBarry Smith if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix"); 5932593348eSBarry Smith 5942593348eSBarry Smith *B = 0; 595b6490206SBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm); 5962593348eSBarry Smith PLogObjectCreate(C); 597b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 5982593348eSBarry Smith PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 599b6490206SBarry Smith C->destroy = MatDestroy_SeqBAIJ; 600b6490206SBarry Smith C->view = MatView_SeqBAIJ; 6012593348eSBarry Smith C->factor = A->factor; 6022593348eSBarry Smith c->row = 0; 6032593348eSBarry Smith c->col = 0; 6042593348eSBarry Smith C->assembled = PETSC_TRUE; 6052593348eSBarry Smith 6062593348eSBarry Smith c->m = a->m; 6072593348eSBarry Smith c->n = a->n; 608b6490206SBarry Smith c->bs = a->bs; 609b6490206SBarry Smith c->mbs = a->mbs; 610b6490206SBarry Smith c->nbs = a->nbs; 6112593348eSBarry Smith 612b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 613b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 614b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 6152593348eSBarry Smith c->imax[i] = a->imax[i]; 6162593348eSBarry Smith c->ilen[i] = a->ilen[i]; 6172593348eSBarry Smith } 6182593348eSBarry Smith 6192593348eSBarry Smith /* allocate the matrix space */ 6202593348eSBarry Smith c->singlemalloc = 1; 621*de6a44a3SBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int)); 6222593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 623*de6a44a3SBarry Smith c->j = (int *) (c->a + nz*bs*bs); 624*de6a44a3SBarry Smith c->i = c->j + nz; 625b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 626b6490206SBarry Smith if (mbs > 0) { 627*de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 6282593348eSBarry Smith if (cpvalues == COPY_VALUES) { 629*de6a44a3SBarry Smith PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar)); 6302593348eSBarry Smith } 6312593348eSBarry Smith } 6322593348eSBarry Smith 633b6490206SBarry Smith PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 6342593348eSBarry Smith c->sorted = a->sorted; 6352593348eSBarry Smith c->roworiented = a->roworiented; 6362593348eSBarry Smith c->nonew = a->nonew; 6372593348eSBarry Smith 6382593348eSBarry Smith if (a->diag) { 639b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 640b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 641b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 6422593348eSBarry Smith c->diag[i] = a->diag[i]; 6432593348eSBarry Smith } 6442593348eSBarry Smith } 6452593348eSBarry Smith else c->diag = 0; 6462593348eSBarry Smith c->nz = a->nz; 6472593348eSBarry Smith c->maxnz = a->maxnz; 6482593348eSBarry Smith c->solve_work = 0; 6492593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 6502593348eSBarry Smith 6512593348eSBarry Smith *B = C; 6522593348eSBarry Smith return 0; 6532593348eSBarry Smith } 6542593348eSBarry Smith 655b6490206SBarry Smith int MatLoad_SeqBAIJ(Viewer bview,MatType type,Mat *A) 6562593348eSBarry Smith { 657b6490206SBarry Smith Mat_SeqBAIJ *a; 6582593348eSBarry Smith Mat B; 659*de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 660b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 66135aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 662*de6a44a3SBarry Smith int *masked, nmask,tmp,ishift,bs2; 663b6490206SBarry Smith Scalar *aa; 6642593348eSBarry Smith PetscObject vobj = (PetscObject) bview; 6652593348eSBarry Smith MPI_Comm comm = vobj->comm; 6662593348eSBarry Smith 66735aab85fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 668*de6a44a3SBarry Smith bs2 = bs*bs; 669b6490206SBarry Smith 6702593348eSBarry Smith MPI_Comm_size(comm,&size); 671b6490206SBarry Smith if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor"); 6722593348eSBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 6732593348eSBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 674*de6a44a3SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object"); 6752593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 6762593348eSBarry Smith 677b6490206SBarry Smith if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices"); 67835aab85fSBarry Smith 67935aab85fSBarry Smith /* 68035aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 68135aab85fSBarry Smith divisible by the blocksize 68235aab85fSBarry Smith */ 683b6490206SBarry Smith mbs = M/bs; 68435aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 68535aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 68635aab85fSBarry Smith else mbs++; 68735aab85fSBarry Smith if (extra_rows) { 68835aab85fSBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize"); 68935aab85fSBarry Smith } 690b6490206SBarry Smith 6912593348eSBarry Smith /* read in row lengths */ 69235aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 6932593348eSBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 69435aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 6952593348eSBarry Smith 696b6490206SBarry Smith /* read in column indices */ 69735aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 698b6490206SBarry Smith ierr = SYRead(fd,jj,nz,SYINT); CHKERRQ(ierr); 69935aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 700b6490206SBarry Smith 701b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 702b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 703b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 70435aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 70535aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 70635aab85fSBarry Smith masked = mask + mbs; 707b6490206SBarry Smith rowcount = 0; nzcount = 0; 708b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 70935aab85fSBarry Smith nmask = 0; 710b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 711b6490206SBarry Smith kmax = rowlengths[rowcount]; 712b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 71335aab85fSBarry Smith tmp = jj[nzcount++]/bs; 71435aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 715b6490206SBarry Smith } 716b6490206SBarry Smith rowcount++; 717b6490206SBarry Smith } 71835aab85fSBarry Smith browlengths[i] += nmask; 71935aab85fSBarry Smith /* zero out the mask elements we set */ 72035aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 721b6490206SBarry Smith } 722b6490206SBarry Smith 7232593348eSBarry Smith /* create our matrix */ 72435aab85fSBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A); 72535aab85fSBarry Smith CHKERRQ(ierr); 7262593348eSBarry Smith B = *A; 727b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 7282593348eSBarry Smith 7292593348eSBarry Smith /* set matrix "i" values */ 730*de6a44a3SBarry Smith a->i[0] = 0; 731b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 732b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 733b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 7342593348eSBarry Smith } 735b6490206SBarry Smith a->nz = 0; 736b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 7372593348eSBarry Smith 738b6490206SBarry Smith /* read in nonzero values */ 73935aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 740b6490206SBarry Smith ierr = SYRead(fd,aa,nz,SYSCALAR); CHKERRQ(ierr); 74135aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 742b6490206SBarry Smith 743b6490206SBarry Smith /* set "a" and "j" values into matrix */ 744b6490206SBarry Smith nzcount = 0; jcount = 0; 745b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 746b6490206SBarry Smith nzcountb = nzcount; 74735aab85fSBarry Smith nmask = 0; 748b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 749b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 750b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 75135aab85fSBarry Smith tmp = jj[nzcount++]/bs; 75235aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 753b6490206SBarry Smith } 754b6490206SBarry Smith rowcount++; 755b6490206SBarry Smith } 756*de6a44a3SBarry Smith /* sort the masked values */ 757*de6a44a3SBarry Smith SYIsort(nmask,masked); 758*de6a44a3SBarry Smith 759b6490206SBarry Smith /* set "j" values into matrix */ 760b6490206SBarry Smith maskcount = 1; 76135aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 76235aab85fSBarry Smith a->j[jcount++] = masked[j]; 763*de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 764b6490206SBarry Smith } 765b6490206SBarry Smith /* set "a" values into matrix */ 766*de6a44a3SBarry Smith ishift = bs2*a->i[i]; 767b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 768b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 769b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 770*de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 771*de6a44a3SBarry Smith block = mask[tmp] - 1; 772*de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 773*de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 774b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 775b6490206SBarry Smith } 776b6490206SBarry Smith } 77735aab85fSBarry Smith /* zero out the mask elements we set */ 77835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 779b6490206SBarry Smith } 78035aab85fSBarry Smith if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix"); 781b6490206SBarry Smith 782b6490206SBarry Smith PetscFree(rowlengths); 783b6490206SBarry Smith PetscFree(browlengths); 784b6490206SBarry Smith PetscFree(aa); 785b6490206SBarry Smith PetscFree(jj); 786b6490206SBarry Smith PetscFree(mask); 787b6490206SBarry Smith 788b6490206SBarry Smith B->assembled = PETSC_TRUE; 789*de6a44a3SBarry Smith 790*de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 791*de6a44a3SBarry Smith if (flg) { 792*de6a44a3SBarry Smith Viewer viewer; 793*de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 794*de6a44a3SBarry Smith ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO,0);CHKERRQ(ierr); 795*de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 796*de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 797*de6a44a3SBarry Smith } 798*de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 799*de6a44a3SBarry Smith if (flg) { 800*de6a44a3SBarry Smith Viewer viewer; 801*de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 802*de6a44a3SBarry Smith ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr); 803*de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 804*de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 805*de6a44a3SBarry Smith } 806*de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr); 807*de6a44a3SBarry Smith if (flg) { 808*de6a44a3SBarry Smith Viewer viewer; 809*de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 810*de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 811*de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 812*de6a44a3SBarry Smith } 813*de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr); 814*de6a44a3SBarry Smith if (flg) { 815*de6a44a3SBarry Smith Viewer viewer; 816*de6a44a3SBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr); 817*de6a44a3SBarry Smith ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_MATLAB,"M");CHKERRQ(ierr); 818*de6a44a3SBarry Smith ierr = MatView(B,viewer); CHKERRQ(ierr); 819*de6a44a3SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 820*de6a44a3SBarry Smith } 821*de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 822*de6a44a3SBarry Smith if (flg) { 823*de6a44a3SBarry Smith Draw win; 824*de6a44a3SBarry Smith ierr = DrawOpenX(B->comm,0,0,0,0,300,300,&win); CHKERRQ(ierr); 825*de6a44a3SBarry Smith ierr = MatView(B,(Viewer)win); CHKERRQ(ierr); 826*de6a44a3SBarry Smith ierr = DrawSyncFlush(win); CHKERRQ(ierr); 827*de6a44a3SBarry Smith ierr = DrawDestroy(win); CHKERRQ(ierr); 828*de6a44a3SBarry Smith } 8292593348eSBarry Smith return 0; 8302593348eSBarry Smith } 8312593348eSBarry Smith 8322593348eSBarry Smith 8332593348eSBarry Smith 834