1b6490206SBarry Smith 22593348eSBarry Smith #ifndef lint 3*35aab85fSBarry Smith static char vcid[] = "$Id: baij.c,v 1.2 1996/02/13 05:21:59 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; 202593348eSBarry Smith int ierr, *ia, *ja,n,*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 n = a->n; 292593348eSBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 302593348eSBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 312593348eSBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 322593348eSBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 332593348eSBarry Smith PetscFree(idx); 342593348eSBarry Smith ISSetPermutation(*rperm); 352593348eSBarry Smith ISSetPermutation(*cperm); 362593348eSBarry Smith ISSetIdentity(*rperm); 372593348eSBarry Smith ISSetIdentity(*cperm); 382593348eSBarry Smith return 0; 392593348eSBarry Smith } 402593348eSBarry Smith 41b6490206SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,a->indexshift,&ia,&ja);CHKERRQ(ierr); 422593348eSBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 432593348eSBarry Smith 442593348eSBarry Smith PetscFree(ia); PetscFree(ja); 452593348eSBarry Smith return 0; 462593348eSBarry Smith } 472593348eSBarry Smith 482593348eSBarry Smith 492593348eSBarry Smith #include "draw.h" 502593348eSBarry Smith #include "pinclude/pviewer.h" 512593348eSBarry Smith #include "sysio.h" 522593348eSBarry Smith 53b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 542593348eSBarry Smith { 55b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 56b6490206SBarry Smith int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l; 57b6490206SBarry Smith Scalar *aa; 582593348eSBarry Smith 592593348eSBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 602593348eSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 612593348eSBarry Smith col_lens[0] = MAT_COOKIE; 622593348eSBarry Smith col_lens[1] = a->m; 632593348eSBarry Smith col_lens[2] = a->n; 64b6490206SBarry Smith col_lens[3] = a->nz*bs*bs; 652593348eSBarry Smith 662593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 67b6490206SBarry Smith count = 0; 68b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 69b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 70b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 71b6490206SBarry Smith } 722593348eSBarry Smith } 732593348eSBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 742593348eSBarry Smith PetscFree(col_lens); 752593348eSBarry Smith 762593348eSBarry Smith /* store column indices (zero start index) */ 77b6490206SBarry Smith jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj); 78b6490206SBarry Smith count = 0; 79b6490206SBarry Smith if (!a->indexshift) { 80b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 81b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 82b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 83b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 84b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 85b6490206SBarry Smith /* printf(" count %d jj %d row %d\n",count-1,jj[count-1],bs*i+j);*/ 862593348eSBarry Smith } 872593348eSBarry Smith } 88b6490206SBarry Smith } 89b6490206SBarry Smith } 90b6490206SBarry Smith } else { 91b6490206SBarry Smith SETERRQ(1,"Not yet done"); 92b6490206SBarry Smith } 93b6490206SBarry Smith ierr = SYWrite(fd,jj,bs*bs*a->nz,SYINT,0); CHKERRQ(ierr); 94b6490206SBarry Smith PetscFree(jj); 952593348eSBarry Smith 962593348eSBarry Smith /* store nonzero values */ 97b6490206SBarry Smith aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa); 98b6490206SBarry Smith count = 0; 99b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 100b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 101b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 102b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 103b6490206SBarry Smith aa[count++] = a->a[bs*bs*k + l*bs + j]; 104b6490206SBarry Smith } 105b6490206SBarry Smith } 106b6490206SBarry Smith } 107b6490206SBarry Smith } 108b6490206SBarry Smith ierr = SYWrite(fd,aa,bs*bs*a->nz,SYSCALAR,0); CHKERRQ(ierr); 109b6490206SBarry Smith PetscFree(aa); 1102593348eSBarry Smith return 0; 1112593348eSBarry Smith } 1122593348eSBarry Smith 113b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 1142593348eSBarry Smith { 115b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 116b6490206SBarry Smith int ierr, i,j, shift = a->indexshift,format,bs = a->bs,k,l; 1172593348eSBarry Smith FILE *fd; 1182593348eSBarry Smith char *outputname; 1192593348eSBarry Smith 1202593348eSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 1212593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 1222593348eSBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 1232593348eSBarry Smith if (format == FILE_FORMAT_INFO) { 1242593348eSBarry Smith /* no need to print additional information */ ; 1252593348eSBarry Smith } 1262593348eSBarry Smith else if (format == FILE_FORMAT_MATLAB) { 127b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported"); 1282593348eSBarry Smith } 1292593348eSBarry Smith else { 130b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 131b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 132b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 133b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 134b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 135b6490206SBarry Smith fprintf(fd," %d %g",bs*a->j[k]+l-shift,a->a[bs*bs*k + l*bs + j]); 1362593348eSBarry Smith } 1372593348eSBarry Smith } 1382593348eSBarry Smith fprintf(fd,"\n"); 1392593348eSBarry Smith } 1402593348eSBarry Smith } 141b6490206SBarry Smith } 1422593348eSBarry Smith fflush(fd); 1432593348eSBarry Smith return 0; 1442593348eSBarry Smith } 1452593348eSBarry Smith 146b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer) 1472593348eSBarry Smith { 1482593348eSBarry Smith Mat A = (Mat) obj; 1492593348eSBarry Smith PetscObject vobj = (PetscObject) viewer; 1502593348eSBarry Smith 1512593348eSBarry Smith if (!viewer) { 1522593348eSBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 1532593348eSBarry Smith } 1542593348eSBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 1552593348eSBarry Smith if (vobj->type == MATLAB_VIEWER) { 156b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported"); 1572593348eSBarry Smith } 1582593348eSBarry Smith else if (vobj->type == ASCII_FILE_VIEWER||vobj->type == ASCII_FILES_VIEWER){ 159b6490206SBarry Smith return MatView_SeqBAIJ_ASCII(A,viewer); 1602593348eSBarry Smith } 1612593348eSBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 162b6490206SBarry Smith return MatView_SeqBAIJ_Binary(A,viewer); 1632593348eSBarry Smith } 1642593348eSBarry Smith } 1652593348eSBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 1662593348eSBarry Smith if (vobj->type == NULLWINDOW) return 0; 167b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported"); 1682593348eSBarry Smith } 1692593348eSBarry Smith return 0; 1702593348eSBarry Smith } 171b6490206SBarry Smith 172b6490206SBarry Smith 173b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A) 1742593348eSBarry Smith { 175b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 176b6490206SBarry Smith PetscMemzero(a->a,a->bs*a->bs*(a->i[a->mbs]+a->indexshift)*sizeof(Scalar)); 1772593348eSBarry Smith return 0; 1782593348eSBarry Smith } 1792593348eSBarry Smith 180b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj) 1812593348eSBarry Smith { 1822593348eSBarry Smith Mat A = (Mat) obj; 183b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1842593348eSBarry Smith 1852593348eSBarry Smith #if defined(PETSC_LOG) 1862593348eSBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1872593348eSBarry Smith #endif 1882593348eSBarry Smith PetscFree(a->a); 1892593348eSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 1902593348eSBarry Smith if (a->diag) PetscFree(a->diag); 1912593348eSBarry Smith if (a->ilen) PetscFree(a->ilen); 1922593348eSBarry Smith if (a->imax) PetscFree(a->imax); 1932593348eSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 1942593348eSBarry Smith PetscFree(a); 1952593348eSBarry Smith PLogObjectDestroy(A); 1962593348eSBarry Smith PetscHeaderDestroy(A); 1972593348eSBarry Smith return 0; 1982593348eSBarry Smith } 1992593348eSBarry Smith 200b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op) 2012593348eSBarry Smith { 202b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2032593348eSBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 2042593348eSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 2052593348eSBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 2062593348eSBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 2072593348eSBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 2082593348eSBarry Smith else if (op == ROWS_SORTED || 2092593348eSBarry Smith op == SYMMETRIC_MATRIX || 2102593348eSBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 2112593348eSBarry Smith op == YES_NEW_DIAGONALS) 212b6490206SBarry Smith PLogInfo((PetscObject)A,"Info:MatSetOption_SeqBAIJ:Option ignored\n"); 2132593348eSBarry Smith else if (op == NO_NEW_DIAGONALS) 214b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");} 2152593348eSBarry Smith else 216b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");} 2172593348eSBarry Smith return 0; 2182593348eSBarry Smith } 2192593348eSBarry Smith 2202593348eSBarry Smith 2212593348eSBarry Smith /* -------------------------------------------------------*/ 2222593348eSBarry Smith /* Should check that shapes of vectors and matrices match */ 2232593348eSBarry Smith /* -------------------------------------------------------*/ 224b6490206SBarry Smith #include "pinclude/plapack.h" 225b6490206SBarry Smith 226b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy) 2272593348eSBarry Smith { 228b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 229b6490206SBarry Smith Scalar *xg,*yg; 230b6490206SBarry Smith register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5; 231b6490206SBarry Smith register Scalar x1,x2,x3,x4,x5; 232b6490206SBarry Smith int mbs = a->mbs, m = a->m, i, *idx,*ii; 233b6490206SBarry Smith int bs = a->bs,j,n,bs2 = bs*bs; 2342593348eSBarry Smith 235b6490206SBarry Smith if (a->indexshift) SETERRQ(PETSC_ERR_SUP,"MatMult_SeqBAIJ:index shift no"); 2362593348eSBarry Smith 237b6490206SBarry Smith VecGetArray(xx,&xg); x = xg; VecGetArray(yy,&yg); y = yg; 238b6490206SBarry Smith PetscMemzero(y,m*sizeof(Scalar)); 239b6490206SBarry Smith x = x; 2402593348eSBarry Smith idx = a->j; 2412593348eSBarry Smith v = a->a; 2422593348eSBarry Smith ii = a->i; 243b6490206SBarry Smith 244b6490206SBarry Smith switch (bs) { 245b6490206SBarry Smith case 1: 2462593348eSBarry Smith for ( i=0; i<m; i++ ) { 2472593348eSBarry Smith n = ii[1] - ii[0]; ii++; 2482593348eSBarry Smith sum = 0.0; 2492593348eSBarry Smith while (n--) sum += *v++ * x[*idx++]; 2502593348eSBarry Smith y[i] = sum; 2512593348eSBarry Smith } 2522593348eSBarry Smith break; 253b6490206SBarry Smith case 2: 254b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 255b6490206SBarry Smith n = ii[1] - ii[0]; ii++; /* number of blocks in row */ 256b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; 257b6490206SBarry Smith for ( j=0; j<n; j++ ) { 258b6490206SBarry Smith xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 259b6490206SBarry Smith sum1 += v[0]*x1 + v[2]*x2; 260b6490206SBarry Smith sum2 += v[1]*x1 + v[3]*x2; 261b6490206SBarry Smith v += 4; 262b6490206SBarry Smith } 263b6490206SBarry Smith y[0] += sum1; y[1] += sum2; 264b6490206SBarry Smith y += 2; 265b6490206SBarry Smith } 266b6490206SBarry Smith break; 267b6490206SBarry Smith case 3: 268b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 269b6490206SBarry Smith n = ii[1] - ii[0]; ii++; /* number of blocks in row */ 270b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 271b6490206SBarry Smith for ( j=0; j<n; j++ ) { 272b6490206SBarry Smith xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 273b6490206SBarry Smith sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 274b6490206SBarry Smith sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 275b6490206SBarry Smith sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 276b6490206SBarry Smith v += 9; 277b6490206SBarry Smith } 278b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; 279b6490206SBarry Smith y += 3; 280b6490206SBarry Smith } 281b6490206SBarry Smith break; 282b6490206SBarry Smith case 4: 283b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 284b6490206SBarry Smith n = ii[1] - ii[0]; ii++; /* number of blocks in row */ 285b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 286b6490206SBarry Smith for ( j=0; j<n; j++ ) { 287b6490206SBarry Smith xb = x + 4*(*idx++); 288b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 289b6490206SBarry Smith sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 290b6490206SBarry Smith sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 291b6490206SBarry Smith sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 292b6490206SBarry Smith sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 293b6490206SBarry Smith v += 16; 294b6490206SBarry Smith } 295b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; 296b6490206SBarry Smith y += 4; 297b6490206SBarry Smith } 298b6490206SBarry Smith break; 299b6490206SBarry Smith case 5: 300b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 301b6490206SBarry Smith n = ii[1] - ii[0]; ii++; /* number of blocks in row */ 302b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 303b6490206SBarry Smith for ( j=0; j<n; j++ ) { 304b6490206SBarry Smith xb = x + 5*(*idx++); 305b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 306b6490206SBarry Smith sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 307b6490206SBarry Smith sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 308b6490206SBarry Smith sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 309b6490206SBarry Smith sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 310b6490206SBarry Smith sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 311b6490206SBarry Smith v += 25; 312b6490206SBarry Smith } 313b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5; 314b6490206SBarry Smith y += 5; 315b6490206SBarry Smith } 316b6490206SBarry Smith break; 317b6490206SBarry Smith /* block sizes larger then 5 by 5 are handled by BLAS */ 318b6490206SBarry Smith default: { 319b6490206SBarry Smith int _One = 1; Scalar _DOne = 1.0; 320b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 321b6490206SBarry Smith n = ii[1] - ii[0]; ii++; /* number of blocks in row */ 322b6490206SBarry Smith for ( j=0; j<n; j++ ) { 323b6490206SBarry Smith xb = x + bs*(*idx++); 324b6490206SBarry Smith /* do dense mat-vec product */ 325b6490206SBarry Smith LAgemv_("N",&bs,&bs,&_DOne,v,&bs,xb,&_One,&_DOne,y,&_One); 326b6490206SBarry Smith v += bs2; 327b6490206SBarry Smith } 328b6490206SBarry Smith y += bs; 3292593348eSBarry Smith } 3302593348eSBarry Smith } 3312593348eSBarry Smith } 332b6490206SBarry Smith PLogFlops(2*bs2*a->nz - m); 3332593348eSBarry Smith return 0; 3342593348eSBarry Smith } 3352593348eSBarry Smith 336b6490206SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 3372593348eSBarry Smith { 338b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 339*35aab85fSBarry Smith *nz = a->bs*a->bs*a->nz; 3402593348eSBarry Smith *nzalloc = a->maxnz; 3412593348eSBarry Smith *mem = (int)A->mem; 3422593348eSBarry Smith return 0; 3432593348eSBarry Smith } 3442593348eSBarry Smith 345b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 346b6490206SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ(Mat,Mat*); 347b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 348b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec); 3492593348eSBarry Smith 350b6490206SBarry Smith static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 3512593348eSBarry Smith { 352b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3532593348eSBarry Smith *m = a->m; *n = a->n; 3542593348eSBarry Smith return 0; 3552593348eSBarry Smith } 3562593348eSBarry Smith 357b6490206SBarry Smith static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 3582593348eSBarry Smith { 359b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3602593348eSBarry Smith *m = 0; *n = a->m; 3612593348eSBarry Smith return 0; 3622593348eSBarry Smith } 363b6490206SBarry Smith 364b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 3652593348eSBarry Smith { 366b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3672593348eSBarry Smith Scalar *v = a->a; 3682593348eSBarry Smith double sum = 0.0; 369b6490206SBarry Smith int i; 3702593348eSBarry Smith 3712593348eSBarry Smith if (type == NORM_FROBENIUS) { 3722593348eSBarry Smith for (i=0; i<a->nz; i++ ) { 3732593348eSBarry Smith #if defined(PETSC_COMPLEX) 3742593348eSBarry Smith sum += real(conj(*v)*(*v)); v++; 3752593348eSBarry Smith #else 3762593348eSBarry Smith sum += (*v)*(*v); v++; 3772593348eSBarry Smith #endif 3782593348eSBarry Smith } 3792593348eSBarry Smith *norm = sqrt(sum); 3802593348eSBarry Smith } 3812593348eSBarry Smith else { 382b6490206SBarry Smith SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet"); 3832593348eSBarry Smith } 3842593348eSBarry Smith return 0; 3852593348eSBarry Smith } 3862593348eSBarry Smith 3872593348eSBarry Smith 3882593348eSBarry Smith /* 3892593348eSBarry Smith note: This can only work for identity for row and col. It would 3902593348eSBarry Smith be good to check this and otherwise generate an error. 3912593348eSBarry Smith */ 392b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 3932593348eSBarry Smith { 394b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 3952593348eSBarry Smith int ierr; 3962593348eSBarry Smith Mat outA; 3972593348eSBarry Smith 398b6490206SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported"); 3992593348eSBarry Smith 4002593348eSBarry Smith outA = inA; 4012593348eSBarry Smith inA->factor = FACTOR_LU; 4022593348eSBarry Smith a->row = row; 4032593348eSBarry Smith a->col = col; 4042593348eSBarry Smith 4052593348eSBarry Smith a->solve_work = (Scalar *) PetscMalloc((a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work); 4062593348eSBarry Smith 407b6490206SBarry Smith /* 4082593348eSBarry Smith if (!a->diag) { 409b6490206SBarry Smith ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 4102593348eSBarry Smith } 411b6490206SBarry Smith ierr = MatLUFactorNumeric_SeqBAIJ(inA,&outA); CHKERRQ(ierr); 412b6490206SBarry Smith */ 4132593348eSBarry Smith return 0; 4142593348eSBarry Smith } 4152593348eSBarry Smith 4162593348eSBarry Smith #include "pinclude/plapack.h" 417b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 4182593348eSBarry Smith { 419b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 420b6490206SBarry Smith int one = 1, totalnz = a->bs*a->bs*a->nz; 421b6490206SBarry Smith BLscal_( &totalnz, alpha, a->a, &one ); 422b6490206SBarry Smith PLogFlops(totalnz); 4232593348eSBarry Smith return 0; 4242593348eSBarry Smith } 4252593348eSBarry Smith 426b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A) 4272593348eSBarry Smith { 4282593348eSBarry Smith static int called = 0; 4292593348eSBarry Smith MPI_Comm comm = A->comm; 4302593348eSBarry Smith 4312593348eSBarry Smith if (called) return 0; else called = 1; 432b6490206SBarry Smith MPIU_printf(comm," Options for MATSeqBAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 4332593348eSBarry Smith MPIU_printf(comm," -mat_lu_pivotthreshold <threshold>\n"); 4342593348eSBarry Smith MPIU_printf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 4352593348eSBarry Smith return 0; 4362593348eSBarry Smith } 4372593348eSBarry Smith /* -------------------------------------------------------------------*/ 438b6490206SBarry Smith static struct _MatOps MatOps = {0, 439b6490206SBarry Smith 0,0, 440b6490206SBarry Smith MatMult_SeqBAIJ,0, 441b6490206SBarry Smith 0,0, 442b6490206SBarry Smith /*MatSolve_SeqBAIJ*/ 0,0, 443b6490206SBarry Smith 0,0, 444b6490206SBarry Smith /*MatLUFactor_SeqBAIJ*/0,0, 445b6490206SBarry Smith 0, 446b6490206SBarry Smith 0, 447b6490206SBarry Smith MatGetInfo_SeqBAIJ,0, 448b6490206SBarry Smith 0,0,MatNorm_SeqBAIJ, 449b6490206SBarry Smith 0,0, 450b6490206SBarry Smith 0, 451b6490206SBarry Smith MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0, 452b6490206SBarry Smith MatGetReordering_SeqBAIJ, 453b6490206SBarry Smith /*MatLUFactorSymbolic_SeqBAIJ */0,/* MatLUFactorNumeric_SeqBAIJ*/ 0,0,0, 454b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 455b6490206SBarry Smith /* MatILUFactorSymbolic_SeqBAIJ */ 0,0, 456b6490206SBarry Smith 0,0,/*MatConvert_SeqBAIJ*/ 0, 457b6490206SBarry Smith 0,0, 458b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 459b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 460b6490206SBarry Smith 0,0, 461b6490206SBarry Smith 0,0, 462b6490206SBarry Smith MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 463b6490206SBarry Smith 0}; 4642593348eSBarry Smith 4652593348eSBarry Smith /*@C 466b6490206SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format 4672593348eSBarry Smith (the default parallel PETSc format). For good matrix assembly performance 4682593348eSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 4692593348eSBarry Smith (or nzz). By setting these parameters accurately, performance can be 4702593348eSBarry Smith increased by more than a factor of 50. 4712593348eSBarry Smith 4722593348eSBarry Smith Input Parameters: 4732593348eSBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 474b6490206SBarry Smith . bs - size of block 4752593348eSBarry Smith . m - number of rows 4762593348eSBarry Smith . n - number of columns 477b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 478b6490206SBarry Smith . nzz - number of block nonzeros per block row or PETSC_NULL 479b6490206SBarry Smith (possibly different for each row) 4802593348eSBarry Smith 4812593348eSBarry Smith Output Parameter: 4822593348eSBarry Smith . A - the matrix 4832593348eSBarry Smith 4842593348eSBarry Smith Notes: 485b6490206SBarry Smith The BAIJ format, is fully compatible with standard Fortran 77 4862593348eSBarry Smith storage. That is, the stored row and column indices can begin at 4872593348eSBarry Smith either one (as in Fortran) or zero. See the users manual for details. 4882593348eSBarry Smith 4892593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 4902593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 4912593348eSBarry Smith allocation. For additional details, see the users manual chapter on 4922593348eSBarry Smith matrices and the file $(PETSC_DIR)/Performance. 4932593348eSBarry Smith 4942593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 4952593348eSBarry Smith @*/ 496b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 4972593348eSBarry Smith { 4982593348eSBarry Smith Mat B; 499b6490206SBarry Smith Mat_SeqBAIJ *b; 500b6490206SBarry Smith int i,len,ierr, flg,mbs = m/bs; 501b6490206SBarry Smith 502b6490206SBarry Smith if (mbs*bs != m) 503b6490206SBarry Smith SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize"); 5042593348eSBarry Smith 5052593348eSBarry Smith *A = 0; 506b6490206SBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm); 5072593348eSBarry Smith PLogObjectCreate(B); 508b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 5092593348eSBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 510b6490206SBarry Smith B->destroy = MatDestroy_SeqBAIJ; 511b6490206SBarry Smith B->view = MatView_SeqBAIJ; 5122593348eSBarry Smith B->factor = 0; 5132593348eSBarry Smith B->lupivotthreshold = 1.0; 5142593348eSBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, \ 5152593348eSBarry Smith &flg); CHKERRQ(ierr); 5162593348eSBarry Smith b->row = 0; 5172593348eSBarry Smith b->col = 0; 5182593348eSBarry Smith b->indexshift = 0; 5192593348eSBarry Smith b->reallocs = 0; 5202593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 5212593348eSBarry Smith if (flg) b->indexshift = -1; 5222593348eSBarry Smith 5232593348eSBarry Smith b->m = m; 524b6490206SBarry Smith b->mbs = mbs; 5252593348eSBarry Smith b->n = n; 526b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 5272593348eSBarry Smith if (nnz == PETSC_NULL) { 5282593348eSBarry Smith if (nz == PETSC_DEFAULT) nz = 10; 5292593348eSBarry Smith else if (nz <= 0) nz = 1; 530b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 531b6490206SBarry Smith nz = nz*mbs; 5322593348eSBarry Smith } 5332593348eSBarry Smith else { 5342593348eSBarry Smith nz = 0; 535b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 5362593348eSBarry Smith } 5372593348eSBarry Smith 5382593348eSBarry Smith /* allocate the matrix space */ 539b6490206SBarry Smith len = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int); 5402593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 541b6490206SBarry Smith PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar)); 542b6490206SBarry Smith b->j = (int *) (b->a + nz*bs*bs); 5432593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 5442593348eSBarry Smith b->i = b->j + nz; 5452593348eSBarry Smith b->singlemalloc = 1; 5462593348eSBarry Smith 5472593348eSBarry Smith b->i[0] = -b->indexshift; 548b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 5492593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 5502593348eSBarry Smith } 5512593348eSBarry Smith 552b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 553b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 554b6490206SBarry Smith PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 555b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 5562593348eSBarry Smith 557b6490206SBarry Smith b->bs = bs; 558b6490206SBarry Smith b->mbs = mbs; 5592593348eSBarry Smith b->nz = 0; 5602593348eSBarry Smith b->maxnz = nz; 5612593348eSBarry Smith b->sorted = 0; 5622593348eSBarry Smith b->roworiented = 1; 5632593348eSBarry Smith b->nonew = 0; 5642593348eSBarry Smith b->diag = 0; 5652593348eSBarry Smith b->solve_work = 0; 5662593348eSBarry Smith b->spptr = 0; 5672593348eSBarry Smith 5682593348eSBarry Smith *A = B; 5692593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 5702593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 5712593348eSBarry Smith return 0; 5722593348eSBarry Smith } 5732593348eSBarry Smith 574b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 5752593348eSBarry Smith { 5762593348eSBarry Smith Mat C; 577b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 578b6490206SBarry Smith int i,len, shift = a->indexshift, mbs = a->mbs, bs = a->bs; 5792593348eSBarry Smith 5802593348eSBarry Smith *B = 0; 581b6490206SBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm); 5822593348eSBarry Smith PLogObjectCreate(C); 583b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 5842593348eSBarry Smith PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 585b6490206SBarry Smith C->destroy = MatDestroy_SeqBAIJ; 586b6490206SBarry Smith C->view = MatView_SeqBAIJ; 5872593348eSBarry Smith C->factor = A->factor; 5882593348eSBarry Smith c->row = 0; 5892593348eSBarry Smith c->col = 0; 5902593348eSBarry Smith c->indexshift = shift; 5912593348eSBarry Smith C->assembled = PETSC_TRUE; 5922593348eSBarry Smith 5932593348eSBarry Smith c->m = a->m; 5942593348eSBarry Smith c->n = a->n; 595b6490206SBarry Smith c->bs = a->bs; 596b6490206SBarry Smith c->mbs = a->mbs; 597b6490206SBarry Smith c->nbs = a->nbs; 5982593348eSBarry Smith 599b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 600b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 601b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 6022593348eSBarry Smith c->imax[i] = a->imax[i]; 6032593348eSBarry Smith c->ilen[i] = a->ilen[i]; 6042593348eSBarry Smith } 6052593348eSBarry Smith 6062593348eSBarry Smith /* allocate the matrix space */ 6072593348eSBarry Smith c->singlemalloc = 1; 608b6490206SBarry Smith len = (mbs+1)*sizeof(int)+(a->i[mbs])*(bs*bs*sizeof(Scalar)+sizeof(int)); 6092593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 610b6490206SBarry Smith c->j = (int *) (c->a + a->i[mbs]*bs*bs + shift); 611b6490206SBarry Smith c->i = c->j + a->i[mbs] + shift; 612b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 613b6490206SBarry Smith if (mbs > 0) { 614b6490206SBarry Smith PetscMemcpy(c->j,a->j,(a->i[mbs]+shift)*sizeof(int)); 6152593348eSBarry Smith if (cpvalues == COPY_VALUES) { 616b6490206SBarry Smith PetscMemcpy(c->a,a->a,(bs*bs*a->i[mbs]+shift)*sizeof(Scalar)); 6172593348eSBarry Smith } 6182593348eSBarry Smith } 6192593348eSBarry Smith 620b6490206SBarry Smith PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 6212593348eSBarry Smith c->sorted = a->sorted; 6222593348eSBarry Smith c->roworiented = a->roworiented; 6232593348eSBarry Smith c->nonew = a->nonew; 6242593348eSBarry Smith 6252593348eSBarry Smith if (a->diag) { 626b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 627b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 628b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 6292593348eSBarry Smith c->diag[i] = a->diag[i]; 6302593348eSBarry Smith } 6312593348eSBarry Smith } 6322593348eSBarry Smith else c->diag = 0; 6332593348eSBarry Smith c->nz = a->nz; 6342593348eSBarry Smith c->maxnz = a->maxnz; 6352593348eSBarry Smith c->solve_work = 0; 6362593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 6372593348eSBarry Smith 6382593348eSBarry Smith *B = C; 6392593348eSBarry Smith return 0; 6402593348eSBarry Smith } 6412593348eSBarry Smith 642b6490206SBarry Smith int MatLoad_SeqBAIJ(Viewer bview,MatType type,Mat *A) 6432593348eSBarry Smith { 644b6490206SBarry Smith Mat_SeqBAIJ *a; 6452593348eSBarry Smith Mat B; 646b6490206SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,shift,bs=1,flg; 647b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 648*35aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 649*35aab85fSBarry Smith int *masked, nmask,tmp; 650b6490206SBarry Smith Scalar *aa; 6512593348eSBarry Smith PetscObject vobj = (PetscObject) bview; 6522593348eSBarry Smith MPI_Comm comm = vobj->comm; 6532593348eSBarry Smith 654*35aab85fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 655b6490206SBarry Smith 6562593348eSBarry Smith MPI_Comm_size(comm,&size); 657b6490206SBarry Smith if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor"); 6582593348eSBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 6592593348eSBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 660b6490206SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object in file"); 6612593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 6622593348eSBarry Smith 663b6490206SBarry Smith if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices"); 664*35aab85fSBarry Smith 665*35aab85fSBarry Smith /* 666*35aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 667*35aab85fSBarry Smith divisible by the blocksize 668*35aab85fSBarry Smith */ 669b6490206SBarry Smith mbs = M/bs; 670*35aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 671*35aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 672*35aab85fSBarry Smith else mbs++; 673*35aab85fSBarry Smith if (extra_rows) { 674*35aab85fSBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize"); 675*35aab85fSBarry Smith } 676b6490206SBarry Smith 6772593348eSBarry Smith /* read in row lengths */ 678*35aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 6792593348eSBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 680*35aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 6812593348eSBarry Smith 682b6490206SBarry Smith /* read in column indices */ 683*35aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 684b6490206SBarry Smith ierr = SYRead(fd,jj,nz,SYINT); CHKERRQ(ierr); 685*35aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 686b6490206SBarry Smith 687b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 688b6490206SBarry Smith browlengths = (int *) PetscMalloc( mbs*sizeof(int) );CHKPTRQ(browlengths); 689b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 690*35aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 691*35aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 692*35aab85fSBarry Smith masked = mask + mbs; 693b6490206SBarry Smith rowcount = 0; nzcount = 0; 694b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 695*35aab85fSBarry Smith nmask = 0; 696b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 697b6490206SBarry Smith kmax = rowlengths[rowcount]; 698b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 699*35aab85fSBarry Smith tmp = jj[nzcount++]/bs; 700*35aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 701b6490206SBarry Smith } 702b6490206SBarry Smith rowcount++; 703b6490206SBarry Smith } 704*35aab85fSBarry Smith browlengths[i] += nmask; 705*35aab85fSBarry Smith /* zero out the mask elements we set */ 706*35aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 707b6490206SBarry Smith } 708b6490206SBarry Smith 7092593348eSBarry Smith /* create our matrix */ 710*35aab85fSBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A); 711*35aab85fSBarry Smith CHKERRQ(ierr); 7122593348eSBarry Smith B = *A; 713b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 7142593348eSBarry Smith shift = a->indexshift; 7152593348eSBarry Smith 7162593348eSBarry Smith /* set matrix "i" values */ 7172593348eSBarry Smith a->i[0] = -shift; 718b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 719b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 720b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 7212593348eSBarry Smith } 722b6490206SBarry Smith a->nz = 0; 723b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 7242593348eSBarry Smith 725b6490206SBarry Smith /* read in nonzero values */ 726*35aab85fSBarry Smith aa = (Scalar *) PetscMalloc( (nz+extra_rows)*sizeof(Scalar) ); CHKPTRQ(aa); 727b6490206SBarry Smith ierr = SYRead(fd,aa,nz,SYSCALAR); CHKERRQ(ierr); 728*35aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 729b6490206SBarry Smith 730b6490206SBarry Smith /* set "a" and "j" values into matrix */ 731b6490206SBarry Smith nzcount = 0; jcount = 0; 732b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 733b6490206SBarry Smith nzcountb = nzcount; 734*35aab85fSBarry Smith nmask = 0; 735b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 736b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 737b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 738*35aab85fSBarry Smith tmp = jj[nzcount++]/bs; 739*35aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 740b6490206SBarry Smith } 741b6490206SBarry Smith rowcount++; 742b6490206SBarry Smith } 743b6490206SBarry Smith /* set "j" values into matrix */ 744b6490206SBarry Smith maskcount = 1; 745*35aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 746*35aab85fSBarry Smith a->j[jcount++] = masked[j]; 747*35aab85fSBarry Smith mask[masked[j]] = maskcount++; /* what nonzero block in this row is j */ 748b6490206SBarry Smith } 749b6490206SBarry Smith /* set "a" values into matrix */ 750b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 751b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 752b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 753b6490206SBarry Smith block = mask[jj[nzcountb]/bs] - 1; 754b6490206SBarry Smith point = jj[nzcountb] - bs*(jj[nzcountb]/bs); 755b6490206SBarry Smith idx = bs*bs*(a->i[i] + block) + j + bs*(point); 756b6490206SBarry Smith /* printf("block row %d subrow %d cold %d block %d idx %d val %g a->i[i] %d point %d\n",i,j,k,block,idx, 757b6490206SBarry Smith aa[nzcountb],a->i[i],point); */ 758b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 759b6490206SBarry Smith } 760b6490206SBarry Smith } 761*35aab85fSBarry Smith /* zero out the mask elements we set */ 762*35aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 763b6490206SBarry Smith } 764*35aab85fSBarry Smith if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix"); 765b6490206SBarry Smith 766b6490206SBarry Smith PetscFree(rowlengths); 767b6490206SBarry Smith PetscFree(browlengths); 768b6490206SBarry Smith PetscFree(aa); 769b6490206SBarry Smith PetscFree(jj); 770b6490206SBarry Smith PetscFree(mask); 771b6490206SBarry Smith 772b6490206SBarry Smith B->assembled = PETSC_TRUE; 7732593348eSBarry Smith return 0; 7742593348eSBarry Smith } 7752593348eSBarry Smith 7762593348eSBarry Smith 7772593348eSBarry Smith 778