1b6490206SBarry Smith 22593348eSBarry Smith #ifndef lint 3*7e67e3f9SSatish Balay static char vcid[] = "$Id: baij.c,v 1.26 1996/04/04 04:07:02 balay Exp balay $"; 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 15bcd2baecSBarry 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; 20bcd2baecSBarry 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 40bcd2baecSBarry Smith MatReorderingRegisterAll(); 41bcd2baecSBarry Smith ishift = a->indexshift; 42bcd2baecSBarry Smith oshift = -MatReorderingIndexShift[(int)type]; 43bcd2baecSBarry Smith if (MatReorderingRequiresSymmetric[(int)type]) { 44bcd2baecSBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja); 45bcd2baecSBarry Smith CHKERRQ(ierr); 46bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 472593348eSBarry Smith PetscFree(ia); PetscFree(ja); 48bcd2baecSBarry Smith } else { 49bcd2baecSBarry Smith if (ishift == oshift) { 50bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 51bcd2baecSBarry Smith } 52bcd2baecSBarry Smith else if (ishift == -1) { 53bcd2baecSBarry Smith /* temporarily subtract 1 from i and j indices */ 54bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 55bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 56bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 57bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 58bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 59bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 60bcd2baecSBarry Smith } else { 61bcd2baecSBarry Smith /* temporarily add 1 to i and j indices */ 62bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 63bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 64bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 65bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 66bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 67bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 68bcd2baecSBarry Smith } 69bcd2baecSBarry 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" 9877c4ece6SBarry Smith #include "sys.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; 103*7e67e3f9SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=bs*bs; 104b6490206SBarry Smith Scalar *aa; 1052593348eSBarry Smith 10690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(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; 111*7e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 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 } 12077c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 1212593348eSBarry Smith PetscFree(col_lens); 1222593348eSBarry Smith 1232593348eSBarry Smith /* store column indices (zero start index) */ 124*7e67e3f9SSatish Balay jj = (int *) PetscMalloc( a->nz*bs2*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 } 135*7e67e3f9SSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,BINARY_INT,0); CHKERRQ(ierr); 136b6490206SBarry Smith PetscFree(jj); 1372593348eSBarry Smith 1382593348eSBarry Smith /* store nonzero values */ 139*7e67e3f9SSatish Balay aa = (Scalar *) PetscMalloc(a->nz*bs2*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++ ) { 145*7e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 146b6490206SBarry Smith } 147b6490206SBarry Smith } 148b6490206SBarry Smith } 149b6490206SBarry Smith } 150*7e67e3f9SSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,BINARY_SCALAR,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; 158*7e67e3f9SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=bs*bs; 1592593348eSBarry Smith FILE *fd; 1602593348eSBarry Smith char *outputname; 1612593348eSBarry Smith 16290ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 1632593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 16490ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 16590ace30eSBarry Smith if (format == ASCII_FORMAT_INFO) { 1662593348eSBarry Smith /* no need to print additional information */ ; 1672593348eSBarry Smith } 16890ace30eSBarry Smith else if (format == ASCII_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) 178*7e67e3f9SSatish Balay if (imag(a->a[bs2*k + l*bs + j]) != 0.0) { 17988685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 180*7e67e3f9SSatish Balay real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 18188685aaeSLois Curfman McInnes } 18288685aaeSLois Curfman McInnes else { 183*7e67e3f9SSatish Balay fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 18488685aaeSLois Curfman McInnes } 18588685aaeSLois Curfman McInnes #else 186*7e67e3f9SSatish Balay fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*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; 20119bcc07fSBarry Smith ViewerType vtype; 20219bcc07fSBarry Smith int ierr; 2032593348eSBarry Smith 2042593348eSBarry Smith if (!viewer) { 20519bcc07fSBarry Smith viewer = STDOUT_VIEWER_SELF; 2062593348eSBarry Smith } 20719bcc07fSBarry Smith 20819bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 20919bcc07fSBarry Smith if (vtype == MATLAB_VIEWER) { 210b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported"); 2112593348eSBarry Smith } 21219bcc07fSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 213b6490206SBarry Smith return MatView_SeqBAIJ_ASCII(A,viewer); 2142593348eSBarry Smith } 21519bcc07fSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 216b6490206SBarry Smith return MatView_SeqBAIJ_Binary(A,viewer); 2172593348eSBarry Smith } 21819bcc07fSBarry Smith else if (vtype == DRAW_VIEWER) { 219b6490206SBarry Smith SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported"); 2202593348eSBarry Smith } 2212593348eSBarry Smith return 0; 2222593348eSBarry Smith } 223b6490206SBarry Smith 224*7e67e3f9SSatish Balay #define CHUNKSIZE 1 225cd0e1443SSatish Balay 226cd0e1443SSatish Balay /* This version has row oriented v */ 227cd0e1443SSatish Balay static int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 228cd0e1443SSatish Balay { 229cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 230cd0e1443SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 231cd0e1443SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 232cd0e1443SSatish Balay int *aj=a->j,nonew=a->nonew,shift=a->indexshift,bs=a->bs,brow,bcol; 233*7e67e3f9SSatish Balay int ridx,cidx,bs2=bs*bs; 234cd0e1443SSatish Balay Scalar *ap,value,*aa=a->a,*bap; 235cd0e1443SSatish Balay 236cd0e1443SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 237cd0e1443SSatish Balay row = im[k]; brow = row/bs; 238cd0e1443SSatish Balay if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row"); 239cd0e1443SSatish Balay if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large"); 240*7e67e3f9SSatish Balay rp = aj + ai[brow] + shift; ap = aa + bs2*ai[brow] + shift; 241cd0e1443SSatish Balay rmax = imax[brow]; nrow = ailen[brow]; 242cd0e1443SSatish Balay low = 0; 243cd0e1443SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 244cd0e1443SSatish Balay if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column"); 245cd0e1443SSatish Balay if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large"); 246cd0e1443SSatish Balay col = in[l] - shift; bcol = col/bs; 247cd0e1443SSatish Balay ridx = row % bs; cidx = col % bs; 248cd0e1443SSatish Balay if (roworiented) { 249cd0e1443SSatish Balay value = *v++; 250cd0e1443SSatish Balay } 251cd0e1443SSatish Balay else { 252cd0e1443SSatish Balay value = v[k + l*m]; 253cd0e1443SSatish Balay } 254cd0e1443SSatish Balay if (!sorted) low = 0; high = nrow; 255cd0e1443SSatish Balay while (high-low > 5) { 256cd0e1443SSatish Balay t = (low+high)/2; 257cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 258cd0e1443SSatish Balay else low = t; 259cd0e1443SSatish Balay } 260cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 261cd0e1443SSatish Balay if (rp[i] > bcol) break; 262cd0e1443SSatish Balay if (rp[i] == bcol) { 263*7e67e3f9SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 264cd0e1443SSatish Balay if (is == ADD_VALUES) *bap += value; 265cd0e1443SSatish Balay else *bap = value; 266cd0e1443SSatish Balay goto noinsert; 267cd0e1443SSatish Balay } 268cd0e1443SSatish Balay } 269cd0e1443SSatish Balay if (nonew) goto noinsert; 270cd0e1443SSatish Balay if (nrow >= rmax) { 271cd0e1443SSatish Balay /* there is no extra room in row, therefore enlarge */ 272cd0e1443SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 273cd0e1443SSatish Balay Scalar *new_a; 274cd0e1443SSatish Balay 275cd0e1443SSatish Balay /* malloc new storage space */ 276*7e67e3f9SSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 277cd0e1443SSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 278*7e67e3f9SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 279cd0e1443SSatish Balay new_i = new_j + new_nz; 280cd0e1443SSatish Balay 281cd0e1443SSatish Balay /* copy over old data into new slots */ 282cd0e1443SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 283*7e67e3f9SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 284*7e67e3f9SSatish Balay PetscMemcpy(new_j,aj,(ai[brow]+nrow+shift)*sizeof(int)); 285cd0e1443SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow - shift); 286cd0e1443SSatish Balay PetscMemcpy(new_j+ai[brow]+shift+nrow+CHUNKSIZE,aj+ai[brow]+shift+nrow, 287cd0e1443SSatish Balay len*sizeof(int)); 288*7e67e3f9SSatish Balay PetscMemcpy(new_a,aa,(ai[brow]+nrow+shift)*bs2*sizeof(Scalar)); 289*7e67e3f9SSatish Balay PetscMemzero(new_a+bs2*(ai[brow]+shift+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 290*7e67e3f9SSatish Balay PetscMemcpy(new_a+bs2*(ai[brow]+shift+nrow+CHUNKSIZE), 291*7e67e3f9SSatish Balay aa+bs2*(ai[brow]+shift+nrow),bs2*len*sizeof(Scalar)); 292cd0e1443SSatish Balay /* free up old matrix storage */ 293cd0e1443SSatish Balay PetscFree(a->a); 294cd0e1443SSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 295cd0e1443SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 296cd0e1443SSatish Balay a->singlemalloc = 1; 297cd0e1443SSatish Balay 298*7e67e3f9SSatish Balay rp = aj + ai[brow] + shift; ap = aa + bs2*ai[brow] + shift; 299cd0e1443SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 300*7e67e3f9SSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 301*7e67e3f9SSatish Balay a->maxnz += bs2*CHUNKSIZE; 302cd0e1443SSatish Balay a->reallocs++; 303*7e67e3f9SSatish Balay a->nz+=bs2; 304cd0e1443SSatish Balay } 305*7e67e3f9SSatish Balay N = nrow++ - 1; 306cd0e1443SSatish Balay /* shift up all the later entries in this row */ 307cd0e1443SSatish Balay for ( ii=N; ii>=i; ii-- ) { 308cd0e1443SSatish Balay rp[ii+1] = rp[ii]; 309*7e67e3f9SSatish Balay PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 310cd0e1443SSatish Balay } 311*7e67e3f9SSatish Balay if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 312cd0e1443SSatish Balay rp[i] = bcol; 313*7e67e3f9SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 314cd0e1443SSatish Balay noinsert:; 315cd0e1443SSatish Balay low = i; 316cd0e1443SSatish Balay } 317cd0e1443SSatish Balay ailen[brow] = nrow; 318cd0e1443SSatish Balay } 319cd0e1443SSatish Balay return 0; 320cd0e1443SSatish Balay } 321cd0e1443SSatish Balay 3220b824a48SSatish Balay static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 3230b824a48SSatish Balay { 3240b824a48SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3250b824a48SSatish Balay *m = a->m; *n = a->n; 3260b824a48SSatish Balay return 0; 3270b824a48SSatish Balay } 3280b824a48SSatish Balay 3290b824a48SSatish Balay static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 3300b824a48SSatish Balay { 3310b824a48SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3320b824a48SSatish Balay *m = 0; *n = a->m; 3330b824a48SSatish Balay return 0; 3340b824a48SSatish Balay } 3350b824a48SSatish Balay 3369867e207SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 3379867e207SSatish Balay { 3389867e207SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 339*7e67e3f9SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 3409867e207SSatish Balay Scalar *aa,*v_i,*aa_i; 3419867e207SSatish Balay 3429867e207SSatish Balay bs = a->bs; 3439867e207SSatish Balay ai = a->i; 3449867e207SSatish Balay aj = a->j; 3459867e207SSatish Balay aa = a->a; 346*7e67e3f9SSatish Balay bs2 = bs*bs; 3479867e207SSatish Balay 3489867e207SSatish Balay if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range"); 3499867e207SSatish Balay 3509867e207SSatish Balay bn = row/bs; /* Block number */ 3519867e207SSatish Balay bp = row % bs; /* Block Position */ 3529867e207SSatish Balay M = ai[bn+1] - ai[bn]; 3539867e207SSatish Balay *nz = bs*M; 3549867e207SSatish Balay 3559867e207SSatish Balay if (v) { 3569867e207SSatish Balay *v = 0; 3579867e207SSatish Balay if (*nz) { 3589867e207SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 3599867e207SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 3609867e207SSatish Balay v_i = *v + i*bs; 361*7e67e3f9SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 362*7e67e3f9SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 3639867e207SSatish Balay } 3649867e207SSatish Balay } 3659867e207SSatish Balay } 3669867e207SSatish Balay 3679867e207SSatish Balay if (idx) { 3689867e207SSatish Balay *idx = 0; 3699867e207SSatish Balay if (*nz) { 3709867e207SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 3719867e207SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 3729867e207SSatish Balay idx_i = *idx + i*bs; 3739867e207SSatish Balay itmp = bs*aj[ai[bn] + i]; 3749867e207SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 3759867e207SSatish Balay } 3769867e207SSatish Balay } 3779867e207SSatish Balay } 3789867e207SSatish Balay return 0; 3799867e207SSatish Balay } 3809867e207SSatish Balay 3819867e207SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 3829867e207SSatish Balay { 3839867e207SSatish Balay if (idx) {if (*idx) PetscFree(*idx);} 3849867e207SSatish Balay if (v) {if (*v) PetscFree(*v);} 3859867e207SSatish Balay return 0; 3869867e207SSatish Balay } 387b6490206SBarry Smith 388584200bdSSatish Balay static int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 389584200bdSSatish Balay { 390584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 391584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 392584200bdSSatish Balay int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 393584200bdSSatish Balay int bs = a->bs, mbs = a->mbs, bs2 = bs*bs; 394584200bdSSatish Balay Scalar *aa = a->a, *ap; 395584200bdSSatish Balay 396584200bdSSatish Balay if (mode == FLUSH_ASSEMBLY) return 0; 397584200bdSSatish Balay 398584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 399584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 400584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 401584200bdSSatish Balay if (fshift) { 402584200bdSSatish Balay ip = aj + ai[i] + shift; ap = aa + bs2*ai[i] + shift; 403584200bdSSatish Balay N = ailen[i]; 404584200bdSSatish Balay for ( j=0; j<N; j++ ) { 405584200bdSSatish Balay ip[j-fshift] = ip[j]; 406*7e67e3f9SSatish Balay PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar)); 407584200bdSSatish Balay } 408584200bdSSatish Balay } 409584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 410584200bdSSatish Balay } 411584200bdSSatish Balay if (mbs) { 412584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 413584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 414584200bdSSatish Balay } 415584200bdSSatish Balay /* reset ilen and imax for each row */ 416584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 417584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 418584200bdSSatish Balay } 419584200bdSSatish Balay a->nz = (ai[m] + shift)*bs2; 420584200bdSSatish Balay 421584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 422584200bdSSatish Balay if (fshift && a->diag) { 423584200bdSSatish Balay PetscFree(a->diag); 424584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 425584200bdSSatish Balay a->diag = 0; 426584200bdSSatish Balay } 427584200bdSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Unneeded storage space %d used %d rows %d\n", 428584200bdSSatish Balay fshift,a->nz,m); 429584200bdSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues %d\n", 430584200bdSSatish Balay a->reallocs); 431584200bdSSatish Balay return 0; 432584200bdSSatish Balay } 433584200bdSSatish Balay 434b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A) 4352593348eSBarry Smith { 436b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 437de6a44a3SBarry Smith PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar)); 4382593348eSBarry Smith return 0; 4392593348eSBarry Smith } 4402593348eSBarry Smith 441b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj) 4422593348eSBarry Smith { 4432593348eSBarry Smith Mat A = (Mat) obj; 444b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4452593348eSBarry Smith 4462593348eSBarry Smith #if defined(PETSC_LOG) 4472593348eSBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 4482593348eSBarry Smith #endif 4492593348eSBarry Smith PetscFree(a->a); 4502593348eSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4512593348eSBarry Smith if (a->diag) PetscFree(a->diag); 4522593348eSBarry Smith if (a->ilen) PetscFree(a->ilen); 4532593348eSBarry Smith if (a->imax) PetscFree(a->imax); 4542593348eSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 455de6a44a3SBarry Smith if (a->mult_work) PetscFree(a->mult_work); 4562593348eSBarry Smith PetscFree(a); 457f2655603SLois Curfman McInnes PLogObjectDestroy(A); 458f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 4592593348eSBarry Smith return 0; 4602593348eSBarry Smith } 4612593348eSBarry Smith 462b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op) 4632593348eSBarry Smith { 464b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4652593348eSBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 4662593348eSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 4672593348eSBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 4682593348eSBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 4692593348eSBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 4702593348eSBarry Smith else if (op == ROWS_SORTED || 4712593348eSBarry Smith op == SYMMETRIC_MATRIX || 4722593348eSBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 4732593348eSBarry Smith op == YES_NEW_DIAGONALS) 47494a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n"); 4752593348eSBarry Smith else if (op == NO_NEW_DIAGONALS) 476b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");} 4772593348eSBarry Smith else 478b6490206SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");} 4792593348eSBarry Smith return 0; 4802593348eSBarry Smith } 4812593348eSBarry Smith 4822593348eSBarry Smith 4832593348eSBarry Smith /* -------------------------------------------------------*/ 4842593348eSBarry Smith /* Should check that shapes of vectors and matrices match */ 4852593348eSBarry Smith /* -------------------------------------------------------*/ 486b6490206SBarry Smith #include "pinclude/plapack.h" 487b6490206SBarry Smith 488b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy) 4892593348eSBarry Smith { 490b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 491b6490206SBarry Smith Scalar *xg,*yg; 492b6490206SBarry Smith register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5; 493b6490206SBarry Smith register Scalar x1,x2,x3,x4,x5; 494b6490206SBarry Smith int mbs = a->mbs, m = a->m, i, *idx,*ii; 495b6490206SBarry Smith int bs = a->bs,j,n,bs2 = bs*bs; 4962593348eSBarry Smith 497b6490206SBarry Smith VecGetArray(xx,&xg); x = xg; VecGetArray(yy,&yg); y = yg; 498b6490206SBarry Smith PetscMemzero(y,m*sizeof(Scalar)); 499b6490206SBarry Smith x = x; 5002593348eSBarry Smith idx = a->j; 5012593348eSBarry Smith v = a->a; 5022593348eSBarry Smith ii = a->i; 503b6490206SBarry Smith 504b6490206SBarry Smith switch (bs) { 505b6490206SBarry Smith case 1: 5062593348eSBarry Smith for ( i=0; i<m; i++ ) { 5072593348eSBarry Smith n = ii[1] - ii[0]; ii++; 5082593348eSBarry Smith sum = 0.0; 5092593348eSBarry Smith while (n--) sum += *v++ * x[*idx++]; 5102593348eSBarry Smith y[i] = sum; 5112593348eSBarry Smith } 5122593348eSBarry Smith break; 513b6490206SBarry Smith case 2: 514b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 515de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 516b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; 517b6490206SBarry Smith for ( j=0; j<n; j++ ) { 518b6490206SBarry Smith xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 519b6490206SBarry Smith sum1 += v[0]*x1 + v[2]*x2; 520b6490206SBarry Smith sum2 += v[1]*x1 + v[3]*x2; 521b6490206SBarry Smith v += 4; 522b6490206SBarry Smith } 523b6490206SBarry Smith y[0] += sum1; y[1] += sum2; 524b6490206SBarry Smith y += 2; 525b6490206SBarry Smith } 526b6490206SBarry Smith break; 527b6490206SBarry Smith case 3: 528b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 529de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 530b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 531b6490206SBarry Smith for ( j=0; j<n; j++ ) { 532b6490206SBarry Smith xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 533b6490206SBarry Smith sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 534b6490206SBarry Smith sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 535b6490206SBarry Smith sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 536b6490206SBarry Smith v += 9; 537b6490206SBarry Smith } 538b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; 539b6490206SBarry Smith y += 3; 540b6490206SBarry Smith } 541b6490206SBarry Smith break; 542b6490206SBarry Smith case 4: 543b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 544de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 545b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 546b6490206SBarry Smith for ( j=0; j<n; j++ ) { 547b6490206SBarry Smith xb = x + 4*(*idx++); 548b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 549b6490206SBarry Smith sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 550b6490206SBarry Smith sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 551b6490206SBarry Smith sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 552b6490206SBarry Smith sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 553b6490206SBarry Smith v += 16; 554b6490206SBarry Smith } 555b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; 556b6490206SBarry Smith y += 4; 557b6490206SBarry Smith } 558b6490206SBarry Smith break; 559b6490206SBarry Smith case 5: 560b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 561de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 562b6490206SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 563b6490206SBarry Smith for ( j=0; j<n; j++ ) { 564b6490206SBarry Smith xb = x + 5*(*idx++); 565b6490206SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 566b6490206SBarry Smith sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 567b6490206SBarry Smith sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 568b6490206SBarry Smith sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 569b6490206SBarry Smith sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 570b6490206SBarry Smith sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 571b6490206SBarry Smith v += 25; 572b6490206SBarry Smith } 573b6490206SBarry Smith y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5; 574b6490206SBarry Smith y += 5; 575b6490206SBarry Smith } 576b6490206SBarry Smith break; 577b6490206SBarry Smith /* block sizes larger then 5 by 5 are handled by BLAS */ 578b6490206SBarry Smith default: { 579de6a44a3SBarry Smith int _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt; 580de6a44a3SBarry Smith if (!a->mult_work) { 581de6a44a3SBarry Smith a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar)); 582de6a44a3SBarry Smith CHKPTRQ(a->mult_work); 583de6a44a3SBarry Smith } 584de6a44a3SBarry Smith work = a->mult_work; 585b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 586de6a44a3SBarry Smith n = ii[1] - ii[0]; ii++; 587de6a44a3SBarry Smith ncols = n*bs; 588de6a44a3SBarry Smith workt = work; 589b6490206SBarry Smith for ( j=0; j<n; j++ ) { 590b6490206SBarry Smith xb = x + bs*(*idx++); 591de6a44a3SBarry Smith for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 592de6a44a3SBarry Smith workt += bs; 593b6490206SBarry Smith } 594de6a44a3SBarry Smith LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One); 595de6a44a3SBarry Smith v += n*bs2; 596b6490206SBarry Smith y += bs; 5972593348eSBarry Smith } 5982593348eSBarry Smith } 5992593348eSBarry Smith } 600b6490206SBarry Smith PLogFlops(2*bs2*a->nz - m); 6012593348eSBarry Smith return 0; 6022593348eSBarry Smith } 6032593348eSBarry Smith 604de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem) 6052593348eSBarry Smith { 606b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 607bcd2baecSBarry Smith if (nz) *nz = a->bs*a->bs*a->nz; 608bcd2baecSBarry Smith if (nza) *nza = a->maxnz; 609bcd2baecSBarry Smith if (mem) *mem = (int)A->mem; 6102593348eSBarry Smith return 0; 6112593348eSBarry Smith } 6122593348eSBarry Smith 61391d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg) 61491d316f6SSatish Balay { 61591d316f6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 61691d316f6SSatish Balay 61791d316f6SSatish Balay if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type"); 61891d316f6SSatish Balay 61991d316f6SSatish Balay /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 62091d316f6SSatish Balay if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| 62191d316f6SSatish Balay (a->nz != b->nz)||(a->indexshift != b->indexshift)) { 62291d316f6SSatish Balay *flg = PETSC_FALSE; return 0; 62391d316f6SSatish Balay } 62491d316f6SSatish Balay 62591d316f6SSatish Balay /* if the a->i are the same */ 62691d316f6SSatish Balay if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) { 62791d316f6SSatish Balay *flg = PETSC_FALSE; return 0; 62891d316f6SSatish Balay } 62991d316f6SSatish Balay 63091d316f6SSatish Balay /* if a->j are the same */ 63191d316f6SSatish Balay if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) { 63291d316f6SSatish Balay *flg = PETSC_FALSE; return 0; 63391d316f6SSatish Balay } 63491d316f6SSatish Balay 63591d316f6SSatish Balay /* if a->a are the same */ 63691d316f6SSatish Balay if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) { 63791d316f6SSatish Balay *flg = PETSC_FALSE; return 0; 63891d316f6SSatish Balay } 63991d316f6SSatish Balay *flg = PETSC_TRUE; 64091d316f6SSatish Balay return 0; 64191d316f6SSatish Balay 64291d316f6SSatish Balay } 64391d316f6SSatish Balay 64491d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 64591d316f6SSatish Balay { 64691d316f6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 647*7e67e3f9SSatish Balay int i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 64817e48fc4SSatish Balay Scalar *x, zero = 0.0,*aa,*aa_j; 64917e48fc4SSatish Balay 65017e48fc4SSatish Balay bs = a->bs; 65117e48fc4SSatish Balay aa = a->a; 65217e48fc4SSatish Balay ai = a->i; 65317e48fc4SSatish Balay aj = a->j; 65417e48fc4SSatish Balay ambs = a->mbs; 655*7e67e3f9SSatish Balay bs2 = bs*bs; 65691d316f6SSatish Balay 65791d316f6SSatish Balay VecSet(&zero,v); 65891d316f6SSatish Balay VecGetArray(v,&x); VecGetLocalSize(v,&n); 6599867e207SSatish Balay if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector"); 66017e48fc4SSatish Balay for ( i=0; i<ambs; i++ ) { 66117e48fc4SSatish Balay for ( j=ai[i]; j<ai[i+1]; j++ ) { 66217e48fc4SSatish Balay if (aj[j] == i) { 66317e48fc4SSatish Balay row = i*bs; 664*7e67e3f9SSatish Balay aa_j = aa+j*bs2; 665*7e67e3f9SSatish Balay for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 66691d316f6SSatish Balay break; 66791d316f6SSatish Balay } 66891d316f6SSatish Balay } 66991d316f6SSatish Balay } 67091d316f6SSatish Balay return 0; 67191d316f6SSatish Balay } 67291d316f6SSatish Balay 6739867e207SSatish Balay static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 6749867e207SSatish Balay { 6759867e207SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6769867e207SSatish Balay Scalar *l,*r,x,*v,*aa,*li,*ri; 677*7e67e3f9SSatish Balay int i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; 6789867e207SSatish Balay 6799867e207SSatish Balay ai = a->i; 6809867e207SSatish Balay aj = a->j; 6819867e207SSatish Balay aa = a->a; 6829867e207SSatish Balay m = a->m; 6839867e207SSatish Balay n = a->n; 6849867e207SSatish Balay bs = a->bs; 6859867e207SSatish Balay mbs = a->mbs; 686*7e67e3f9SSatish Balay bs2 = bs*bs; 6879867e207SSatish Balay if (ll) { 6889867e207SSatish Balay VecGetArray(ll,&l); VecGetSize(ll,&lm); 6899867e207SSatish Balay if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length"); 6909867e207SSatish Balay for ( i=0; i<mbs; i++ ) { /* for each block row */ 6919867e207SSatish Balay M = ai[i+1] - ai[i]; 6929867e207SSatish Balay li = l + i*bs; 693*7e67e3f9SSatish Balay v = aa + bs2*ai[i]; 6949867e207SSatish Balay for ( j=0; j<M; j++ ) { /* for each block */ 695*7e67e3f9SSatish Balay for ( k=0; k<bs2; k++ ) { 6969867e207SSatish Balay (*v++) *= li[k%bs]; 6979867e207SSatish Balay } 6989867e207SSatish Balay } 6999867e207SSatish Balay } 7009867e207SSatish Balay } 7019867e207SSatish Balay 7029867e207SSatish Balay if (rr) { 7039867e207SSatish Balay VecGetArray(rr,&r); VecGetSize(rr,&rn); 7049867e207SSatish Balay if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length"); 7059867e207SSatish Balay for ( i=0; i<mbs; i++ ) { /* for each block row */ 7069867e207SSatish Balay M = ai[i+1] - ai[i]; 707*7e67e3f9SSatish Balay v = aa + bs2*ai[i]; 7089867e207SSatish Balay for ( j=0; j<M; j++ ) { /* for each block */ 7099867e207SSatish Balay ri = r + bs*aj[ai[i]+j]; 7109867e207SSatish Balay for ( k=0; k<bs; k++ ) { 7119867e207SSatish Balay x = ri[k]; 7129867e207SSatish Balay for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x; 7139867e207SSatish Balay } 7149867e207SSatish Balay } 7159867e207SSatish Balay } 7169867e207SSatish Balay } 7179867e207SSatish Balay return 0; 7189867e207SSatish Balay } 7199867e207SSatish Balay 7209867e207SSatish Balay 721b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 722b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 723b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec); 7247fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 7257fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 7267fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 7277fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 7287fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 7297fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 7302593348eSBarry Smith 731b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 7322593348eSBarry Smith { 733b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 7342593348eSBarry Smith Scalar *v = a->a; 7352593348eSBarry Smith double sum = 0.0; 736b6490206SBarry Smith int i; 7372593348eSBarry Smith 7382593348eSBarry Smith if (type == NORM_FROBENIUS) { 7392593348eSBarry Smith for (i=0; i<a->nz; i++ ) { 7402593348eSBarry Smith #if defined(PETSC_COMPLEX) 7412593348eSBarry Smith sum += real(conj(*v)*(*v)); v++; 7422593348eSBarry Smith #else 7432593348eSBarry Smith sum += (*v)*(*v); v++; 7442593348eSBarry Smith #endif 7452593348eSBarry Smith } 7462593348eSBarry Smith *norm = sqrt(sum); 7472593348eSBarry Smith } 7482593348eSBarry Smith else { 749b6490206SBarry Smith SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet"); 7502593348eSBarry Smith } 7512593348eSBarry Smith return 0; 7522593348eSBarry Smith } 7532593348eSBarry Smith 7542593348eSBarry Smith /* 7552593348eSBarry Smith note: This can only work for identity for row and col. It would 7562593348eSBarry Smith be good to check this and otherwise generate an error. 7572593348eSBarry Smith */ 758b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 7592593348eSBarry Smith { 760b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 7612593348eSBarry Smith Mat outA; 762de6a44a3SBarry Smith int ierr; 7632593348eSBarry Smith 764b6490206SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported"); 7652593348eSBarry Smith 7662593348eSBarry Smith outA = inA; 7672593348eSBarry Smith inA->factor = FACTOR_LU; 7682593348eSBarry Smith a->row = row; 7692593348eSBarry Smith a->col = col; 7702593348eSBarry Smith 771de6a44a3SBarry Smith a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 7722593348eSBarry Smith 7732593348eSBarry Smith if (!a->diag) { 774b6490206SBarry Smith ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 7752593348eSBarry Smith } 7767fc0212eSBarry Smith ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 7772593348eSBarry Smith return 0; 7782593348eSBarry Smith } 7792593348eSBarry Smith 780b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 7812593348eSBarry Smith { 782b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 783b6490206SBarry Smith int one = 1, totalnz = a->bs*a->bs*a->nz; 784b6490206SBarry Smith BLscal_( &totalnz, alpha, a->a, &one ); 785b6490206SBarry Smith PLogFlops(totalnz); 7862593348eSBarry Smith return 0; 7872593348eSBarry Smith } 7882593348eSBarry Smith 789*7e67e3f9SSatish Balay static int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 790*7e67e3f9SSatish Balay { 791*7e67e3f9SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 792*7e67e3f9SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 793*7e67e3f9SSatish Balay int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 794*7e67e3f9SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=bs*bs; 795*7e67e3f9SSatish Balay Scalar *ap, *aa = a->a, zero = 0.0; 796*7e67e3f9SSatish Balay 797*7e67e3f9SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 798*7e67e3f9SSatish Balay row = im[k]; brow = row/bs; 799*7e67e3f9SSatish Balay if (row < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative row"); 800*7e67e3f9SSatish Balay if (row >= a->m) SETERRQ(1,"MatGetValues_SeqBAIJ:Row too large"); 801*7e67e3f9SSatish Balay rp = aj + ai[brow] + shift; ap = aa + bs2*ai[brow] + shift; 802*7e67e3f9SSatish Balay nrow = ailen[brow]; 803*7e67e3f9SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 804*7e67e3f9SSatish Balay if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative column"); 805*7e67e3f9SSatish Balay if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqBAIJ:Column too large"); 806*7e67e3f9SSatish Balay col = in[l] - shift; 807*7e67e3f9SSatish Balay bcol = col/bs; 808*7e67e3f9SSatish Balay cidx = col%bs; 809*7e67e3f9SSatish Balay ridx = row%bs; 810*7e67e3f9SSatish Balay high = nrow; 811*7e67e3f9SSatish Balay low = 0; /* assume unsorted */ 812*7e67e3f9SSatish Balay while (high-low > 5) { 813*7e67e3f9SSatish Balay t = (low+high)/2; 814*7e67e3f9SSatish Balay if (rp[t] > bcol) high = t; 815*7e67e3f9SSatish Balay else low = t; 816*7e67e3f9SSatish Balay } 817*7e67e3f9SSatish Balay for ( i=low; i<high; i++ ) { 818*7e67e3f9SSatish Balay if (rp[i] > bcol) break; 819*7e67e3f9SSatish Balay if (rp[i] == bcol) { 820*7e67e3f9SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 821*7e67e3f9SSatish Balay goto finished; 822*7e67e3f9SSatish Balay } 823*7e67e3f9SSatish Balay } 824*7e67e3f9SSatish Balay *v++ = zero; 825*7e67e3f9SSatish Balay finished:; 826*7e67e3f9SSatish Balay } 827*7e67e3f9SSatish Balay } 828*7e67e3f9SSatish Balay return 0; 829*7e67e3f9SSatish Balay } 830*7e67e3f9SSatish Balay 831b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A) 8322593348eSBarry Smith { 8332593348eSBarry Smith static int called = 0; 8342593348eSBarry Smith 8352593348eSBarry Smith if (called) return 0; else called = 1; 8362593348eSBarry Smith return 0; 8372593348eSBarry Smith } 8382593348eSBarry Smith /* -------------------------------------------------------------------*/ 839cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ, 8409867e207SSatish Balay MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ, 841b6490206SBarry Smith MatMult_SeqBAIJ,0, 842b6490206SBarry Smith 0,0, 843de6a44a3SBarry Smith MatSolve_SeqBAIJ,0, 844b6490206SBarry Smith 0,0, 845de6a44a3SBarry Smith MatLUFactor_SeqBAIJ,0, 846b6490206SBarry Smith 0, 847b6490206SBarry Smith 0, 84817e48fc4SSatish Balay MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ, 8499867e207SSatish Balay MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ, 850584200bdSSatish Balay 0,MatAssemblyEnd_SeqBAIJ, 851b6490206SBarry Smith 0, 852b6490206SBarry Smith MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0, 853b6490206SBarry Smith MatGetReordering_SeqBAIJ, 8547fc0212eSBarry Smith MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 855b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 856de6a44a3SBarry Smith MatILUFactorSymbolic_SeqBAIJ,0, 857b6490206SBarry Smith 0,0,/* MatConvert_SeqBAIJ */ 0, 858b6490206SBarry Smith 0,0, 859b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 860b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 861b6490206SBarry Smith 0,0, 862*7e67e3f9SSatish Balay MatGetValues_SeqBAIJ,0, 863b6490206SBarry Smith MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 864b6490206SBarry Smith 0}; 8652593348eSBarry Smith 8662593348eSBarry Smith /*@C 867b6490206SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format 8682593348eSBarry Smith (the default parallel PETSc format). For good matrix assembly performance 8692593348eSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8702593348eSBarry Smith (or nzz). By setting these parameters accurately, performance can be 8712593348eSBarry Smith increased by more than a factor of 50. 8722593348eSBarry Smith 8732593348eSBarry Smith Input Parameters: 8742593348eSBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 875b6490206SBarry Smith . bs - size of block 8762593348eSBarry Smith . m - number of rows 8772593348eSBarry Smith . n - number of columns 878b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 879b6490206SBarry Smith . nzz - number of block nonzeros per block row or PETSC_NULL 880b6490206SBarry Smith (possibly different for each row) 8812593348eSBarry Smith 8822593348eSBarry Smith Output Parameter: 8832593348eSBarry Smith . A - the matrix 8842593348eSBarry Smith 8852593348eSBarry Smith Notes: 886b6490206SBarry Smith The BAIJ format, is fully compatible with standard Fortran 77 8872593348eSBarry Smith storage. That is, the stored row and column indices can begin at 8882593348eSBarry Smith either one (as in Fortran) or zero. See the users manual for details. 8892593348eSBarry Smith 8902593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 8912593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 8922593348eSBarry Smith allocation. For additional details, see the users manual chapter on 8932593348eSBarry Smith matrices and the file $(PETSC_DIR)/Performance. 8942593348eSBarry Smith 8952593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 8962593348eSBarry Smith @*/ 897b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 8982593348eSBarry Smith { 8992593348eSBarry Smith Mat B; 900b6490206SBarry Smith Mat_SeqBAIJ *b; 901*7e67e3f9SSatish Balay int i,len,ierr, flg,mbs = m/bs,bs2; 902b6490206SBarry Smith 903b6490206SBarry Smith if (mbs*bs != m) 904b6490206SBarry Smith SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize"); 9052593348eSBarry Smith 9062593348eSBarry Smith *A = 0; 907*7e67e3f9SSatish Balay bs2 = bs*bs; 908b6490206SBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm); 9092593348eSBarry Smith PLogObjectCreate(B); 910b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 9112593348eSBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 9127fc0212eSBarry Smith switch (bs) { 9137fc0212eSBarry Smith case 1: 9147fc0212eSBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 9157fc0212eSBarry Smith break; 9164eeb42bcSBarry Smith case 2: 9174eeb42bcSBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 9186d84be18SBarry Smith break; 919f327dd97SBarry Smith case 3: 920f327dd97SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 9214eeb42bcSBarry Smith break; 9226d84be18SBarry Smith case 4: 9236d84be18SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 9246d84be18SBarry Smith break; 9256d84be18SBarry Smith case 5: 9266d84be18SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 9276d84be18SBarry Smith break; 9287fc0212eSBarry Smith } 929b6490206SBarry Smith B->destroy = MatDestroy_SeqBAIJ; 930b6490206SBarry Smith B->view = MatView_SeqBAIJ; 9312593348eSBarry Smith B->factor = 0; 9322593348eSBarry Smith B->lupivotthreshold = 1.0; 9332593348eSBarry Smith b->row = 0; 9342593348eSBarry Smith b->col = 0; 9352593348eSBarry Smith b->reallocs = 0; 9362593348eSBarry Smith 9372593348eSBarry Smith b->m = m; 938b6490206SBarry Smith b->mbs = mbs; 9392593348eSBarry Smith b->n = n; 940b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 9412593348eSBarry Smith if (nnz == PETSC_NULL) { 942*7e67e3f9SSatish Balay if (nz == PETSC_DEFAULT) nz = CHUNKSIZE; 9432593348eSBarry Smith else if (nz <= 0) nz = 1; 944b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 945b6490206SBarry Smith nz = nz*mbs; 9462593348eSBarry Smith } 9472593348eSBarry Smith else { 9482593348eSBarry Smith nz = 0; 949b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 9502593348eSBarry Smith } 9512593348eSBarry Smith 9522593348eSBarry Smith /* allocate the matrix space */ 953*7e67e3f9SSatish Balay len = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int); 9542593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 955*7e67e3f9SSatish Balay PetscMemzero(b->a,nz*bs2*sizeof(Scalar)); 956*7e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 9572593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 9582593348eSBarry Smith b->i = b->j + nz; 9592593348eSBarry Smith b->singlemalloc = 1; 9602593348eSBarry Smith 961de6a44a3SBarry Smith b->i[0] = 0; 962b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 9632593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 9642593348eSBarry Smith } 9652593348eSBarry Smith 966b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 967b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 968b6490206SBarry Smith PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 969b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 9702593348eSBarry Smith 971b6490206SBarry Smith b->bs = bs; 972b6490206SBarry Smith b->mbs = mbs; 9732593348eSBarry Smith b->nz = 0; 9742593348eSBarry Smith b->maxnz = nz; 9752593348eSBarry Smith b->sorted = 0; 9762593348eSBarry Smith b->roworiented = 1; 9772593348eSBarry Smith b->nonew = 0; 9782593348eSBarry Smith b->diag = 0; 9792593348eSBarry Smith b->solve_work = 0; 980de6a44a3SBarry Smith b->mult_work = 0; 9812593348eSBarry Smith b->spptr = 0; 9822593348eSBarry Smith 9832593348eSBarry Smith *A = B; 9842593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 9852593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 9862593348eSBarry Smith return 0; 9872593348eSBarry Smith } 9882593348eSBarry Smith 989b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 9902593348eSBarry Smith { 9912593348eSBarry Smith Mat C; 992b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 993*7e67e3f9SSatish Balay int i,len, mbs = a->mbs, bs = a->bs,nz = a->nz,bs2 = bs*bs; 994de6a44a3SBarry Smith 995de6a44a3SBarry Smith if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix"); 9962593348eSBarry Smith 9972593348eSBarry Smith *B = 0; 998b6490206SBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm); 9992593348eSBarry Smith PLogObjectCreate(C); 1000b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 10012593348eSBarry Smith PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1002b6490206SBarry Smith C->destroy = MatDestroy_SeqBAIJ; 1003b6490206SBarry Smith C->view = MatView_SeqBAIJ; 10042593348eSBarry Smith C->factor = A->factor; 10052593348eSBarry Smith c->row = 0; 10062593348eSBarry Smith c->col = 0; 10072593348eSBarry Smith C->assembled = PETSC_TRUE; 10082593348eSBarry Smith 10092593348eSBarry Smith c->m = a->m; 10102593348eSBarry Smith c->n = a->n; 1011b6490206SBarry Smith c->bs = a->bs; 1012b6490206SBarry Smith c->mbs = a->mbs; 1013b6490206SBarry Smith c->nbs = a->nbs; 10142593348eSBarry Smith 1015b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1016b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1017b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 10182593348eSBarry Smith c->imax[i] = a->imax[i]; 10192593348eSBarry Smith c->ilen[i] = a->ilen[i]; 10202593348eSBarry Smith } 10212593348eSBarry Smith 10222593348eSBarry Smith /* allocate the matrix space */ 10232593348eSBarry Smith c->singlemalloc = 1; 1024*7e67e3f9SSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int)); 10252593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1026*7e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1027de6a44a3SBarry Smith c->i = c->j + nz; 1028b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1029b6490206SBarry Smith if (mbs > 0) { 1030de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 10312593348eSBarry Smith if (cpvalues == COPY_VALUES) { 1032*7e67e3f9SSatish Balay PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar)); 10332593348eSBarry Smith } 10342593348eSBarry Smith } 10352593348eSBarry Smith 1036b6490206SBarry Smith PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 10372593348eSBarry Smith c->sorted = a->sorted; 10382593348eSBarry Smith c->roworiented = a->roworiented; 10392593348eSBarry Smith c->nonew = a->nonew; 10402593348eSBarry Smith 10412593348eSBarry Smith if (a->diag) { 1042b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1043b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1044b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 10452593348eSBarry Smith c->diag[i] = a->diag[i]; 10462593348eSBarry Smith } 10472593348eSBarry Smith } 10482593348eSBarry Smith else c->diag = 0; 10492593348eSBarry Smith c->nz = a->nz; 10502593348eSBarry Smith c->maxnz = a->maxnz; 10512593348eSBarry Smith c->solve_work = 0; 10522593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 10537fc0212eSBarry Smith c->mult_work = 0; 10542593348eSBarry Smith *B = C; 10552593348eSBarry Smith return 0; 10562593348eSBarry Smith } 10572593348eSBarry Smith 105819bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 10592593348eSBarry Smith { 1060b6490206SBarry Smith Mat_SeqBAIJ *a; 10612593348eSBarry Smith Mat B; 1062de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1063b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 106435aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1065de6a44a3SBarry Smith int *masked, nmask,tmp,ishift,bs2; 1066b6490206SBarry Smith Scalar *aa; 106719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 10682593348eSBarry Smith 106935aab85fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 1070de6a44a3SBarry Smith bs2 = bs*bs; 1071b6490206SBarry Smith 10722593348eSBarry Smith MPI_Comm_size(comm,&size); 1073b6490206SBarry Smith if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor"); 107490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 107577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 1076de6a44a3SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object"); 10772593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 10782593348eSBarry Smith 1079b6490206SBarry Smith if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices"); 108035aab85fSBarry Smith 108135aab85fSBarry Smith /* 108235aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 108335aab85fSBarry Smith divisible by the blocksize 108435aab85fSBarry Smith */ 1085b6490206SBarry Smith mbs = M/bs; 108635aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 108735aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 108835aab85fSBarry Smith else mbs++; 108935aab85fSBarry Smith if (extra_rows) { 109035aab85fSBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize"); 109135aab85fSBarry Smith } 1092b6490206SBarry Smith 10932593348eSBarry Smith /* read in row lengths */ 109435aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 109577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 109635aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 10972593348eSBarry Smith 1098b6490206SBarry Smith /* read in column indices */ 109935aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 110077c4ece6SBarry Smith ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr); 110135aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1102b6490206SBarry Smith 1103b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1104b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1105b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 110635aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 110735aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 110835aab85fSBarry Smith masked = mask + mbs; 1109b6490206SBarry Smith rowcount = 0; nzcount = 0; 1110b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 111135aab85fSBarry Smith nmask = 0; 1112b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1113b6490206SBarry Smith kmax = rowlengths[rowcount]; 1114b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 111535aab85fSBarry Smith tmp = jj[nzcount++]/bs; 111635aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1117b6490206SBarry Smith } 1118b6490206SBarry Smith rowcount++; 1119b6490206SBarry Smith } 112035aab85fSBarry Smith browlengths[i] += nmask; 112135aab85fSBarry Smith /* zero out the mask elements we set */ 112235aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1123b6490206SBarry Smith } 1124b6490206SBarry Smith 11252593348eSBarry Smith /* create our matrix */ 112635aab85fSBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A); 112735aab85fSBarry Smith CHKERRQ(ierr); 11282593348eSBarry Smith B = *A; 1129b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 11302593348eSBarry Smith 11312593348eSBarry Smith /* set matrix "i" values */ 1132de6a44a3SBarry Smith a->i[0] = 0; 1133b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1134b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1135b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 11362593348eSBarry Smith } 11379867e207SSatish Balay a->indexshift = 0; 1138b6490206SBarry Smith a->nz = 0; 1139b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 11402593348eSBarry Smith 1141b6490206SBarry Smith /* read in nonzero values */ 114235aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 114377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr); 114435aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1145b6490206SBarry Smith 1146b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1147b6490206SBarry Smith nzcount = 0; jcount = 0; 1148b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1149b6490206SBarry Smith nzcountb = nzcount; 115035aab85fSBarry Smith nmask = 0; 1151b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1152b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1153b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 115435aab85fSBarry Smith tmp = jj[nzcount++]/bs; 115535aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1156b6490206SBarry Smith } 1157b6490206SBarry Smith rowcount++; 1158b6490206SBarry Smith } 1159de6a44a3SBarry Smith /* sort the masked values */ 116077c4ece6SBarry Smith PetscSortInt(nmask,masked); 1161de6a44a3SBarry Smith 1162b6490206SBarry Smith /* set "j" values into matrix */ 1163b6490206SBarry Smith maskcount = 1; 116435aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 116535aab85fSBarry Smith a->j[jcount++] = masked[j]; 1166de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1167b6490206SBarry Smith } 1168b6490206SBarry Smith /* set "a" values into matrix */ 1169de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1170b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1171b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1172b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1173de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1174de6a44a3SBarry Smith block = mask[tmp] - 1; 1175de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1176de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1177b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1178b6490206SBarry Smith } 1179b6490206SBarry Smith } 118035aab85fSBarry Smith /* zero out the mask elements we set */ 118135aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1182b6490206SBarry Smith } 118335aab85fSBarry Smith if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix"); 1184b6490206SBarry Smith 1185b6490206SBarry Smith PetscFree(rowlengths); 1186b6490206SBarry Smith PetscFree(browlengths); 1187b6490206SBarry Smith PetscFree(aa); 1188b6490206SBarry Smith PetscFree(jj); 1189b6490206SBarry Smith PetscFree(mask); 1190b6490206SBarry Smith 1191b6490206SBarry Smith B->assembled = PETSC_TRUE; 1192de6a44a3SBarry Smith 1193de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 1194de6a44a3SBarry Smith if (flg) { 119519bcc07fSBarry Smith Viewer tviewer; 119619bcc07fSBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 119790ace30eSBarry Smith ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr); 119819bcc07fSBarry Smith ierr = MatView(B,tviewer); CHKERRQ(ierr); 119919bcc07fSBarry Smith ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1200de6a44a3SBarry Smith } 1201de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 1202de6a44a3SBarry Smith if (flg) { 120319bcc07fSBarry Smith Viewer tviewer; 120419bcc07fSBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 120590ace30eSBarry Smith ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr); 120619bcc07fSBarry Smith ierr = MatView(B,tviewer); CHKERRQ(ierr); 120719bcc07fSBarry Smith ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1208de6a44a3SBarry Smith } 1209de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr); 1210de6a44a3SBarry Smith if (flg) { 121119bcc07fSBarry Smith Viewer tviewer; 121219bcc07fSBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 121319bcc07fSBarry Smith ierr = MatView(B,tviewer); CHKERRQ(ierr); 121419bcc07fSBarry Smith ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1215de6a44a3SBarry Smith } 1216de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr); 1217de6a44a3SBarry Smith if (flg) { 121819bcc07fSBarry Smith Viewer tviewer; 121919bcc07fSBarry Smith ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 122090ace30eSBarry Smith ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr); 122119bcc07fSBarry Smith ierr = MatView(B,tviewer); CHKERRQ(ierr); 122219bcc07fSBarry Smith ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1223de6a44a3SBarry Smith } 1224de6a44a3SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 1225de6a44a3SBarry Smith if (flg) { 122619bcc07fSBarry Smith Viewer tviewer; 122719bcc07fSBarry Smith ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr); 122819bcc07fSBarry Smith ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr); 122919bcc07fSBarry Smith ierr = ViewerFlush(tviewer); CHKERRQ(ierr); 123019bcc07fSBarry Smith ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1231de6a44a3SBarry Smith } 12322593348eSBarry Smith return 0; 12332593348eSBarry Smith } 12342593348eSBarry Smith 12352593348eSBarry Smith 12362593348eSBarry Smith 1237