1*b57c8cebSBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*b57c8cebSBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.226 1997/11/07 01:58:23 bsmith Exp bsmith $"; 4cb512458SBarry Smith #endif 58a729477SBarry Smith 63369ce9aSBarry Smith #include "pinclude/pviewer.h" 770f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 8f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 9d9942c19SSatish Balay #include "src/inline/spops.h" 108a729477SBarry Smith 119e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column 129e25ed09SBarry Smith number to the local number in the off-diagonal part of the local 139e25ed09SBarry Smith storage of the matrix. This is done in a non scable way since the 149e25ed09SBarry Smith length of colmap equals the global matrix length. 159e25ed09SBarry Smith */ 165615d1e5SSatish Balay #undef __FUNC__ 17d4bb536fSBarry Smith #define __FUNC__ "CreateColmap_MPIAIJ_Private" 180a198c4cSBarry Smith int CreateColmap_MPIAIJ_Private(Mat mat) 199e25ed09SBarry Smith { 2044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 21ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data; 22905e6a2fSBarry Smith int n = B->n,i; 23dbb450caSBarry Smith 243a40ed3dSBarry Smith PetscFunctionBegin; 25758f045eSSatish Balay aij->colmap = (int *) PetscMalloc((aij->N+1)*sizeof(int));CHKPTRQ(aij->colmap); 26464493b3SBarry Smith PLogObjectMemory(mat,aij->N*sizeof(int)); 27cddf8d76SBarry Smith PetscMemzero(aij->colmap,aij->N*sizeof(int)); 28905e6a2fSBarry Smith for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1; 293a40ed3dSBarry Smith PetscFunctionReturn(0); 309e25ed09SBarry Smith } 319e25ed09SBarry Smith 322493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat); 332493cbb0SBarry Smith 340520107fSSatish Balay #define CHUNKSIZE 15 3530770e4dSSatish Balay #define MatSetValues_SeqAIJ_A_Private(row,col,value,addv) \ 360520107fSSatish Balay { \ 370520107fSSatish Balay \ 380520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 3930770e4dSSatish Balay rmax = aimax[row]; nrow = ailen[row]; \ 40f5e9677aSSatish Balay col1 = col - shift; \ 41f5e9677aSSatish Balay \ 42ba4e3ef2SSatish Balay low = 0; high = nrow; \ 43ba4e3ef2SSatish Balay while (high-low > 5) { \ 44ba4e3ef2SSatish Balay t = (low+high)/2; \ 45ba4e3ef2SSatish Balay if (rp[t] > col) high = t; \ 46ba4e3ef2SSatish Balay else low = t; \ 47ba4e3ef2SSatish Balay } \ 480520107fSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 49f5e9677aSSatish Balay if (rp[_i] > col1) break; \ 50f5e9677aSSatish Balay if (rp[_i] == col1) { \ 510520107fSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 520520107fSSatish Balay else ap[_i] = value; \ 5330770e4dSSatish Balay goto a_noinsert; \ 540520107fSSatish Balay } \ 550520107fSSatish Balay } \ 5689280ab3SLois Curfman McInnes if (nonew == 1) goto a_noinsert; \ 5711ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 580520107fSSatish Balay if (nrow >= rmax) { \ 590520107fSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 600520107fSSatish Balay int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \ 610520107fSSatish Balay Scalar *new_a; \ 620520107fSSatish Balay \ 6311ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 6489280ab3SLois Curfman McInnes \ 650520107fSSatish Balay /* malloc new storage space */ \ 660520107fSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \ 670520107fSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 680520107fSSatish Balay new_j = (int *) (new_a + new_nz); \ 690520107fSSatish Balay new_i = new_j + new_nz; \ 700520107fSSatish Balay \ 710520107fSSatish Balay /* copy over old data into new slots */ \ 720520107fSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \ 730520107fSSatish Balay for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \ 740520107fSSatish Balay PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \ 750520107fSSatish Balay len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \ 760520107fSSatish Balay PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \ 770520107fSSatish Balay len*sizeof(int)); \ 780520107fSSatish Balay PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \ 790520107fSSatish Balay PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \ 800520107fSSatish Balay len*sizeof(Scalar)); \ 810520107fSSatish Balay /* free up old matrix storage */ \ 82f5e9677aSSatish Balay \ 830520107fSSatish Balay PetscFree(a->a); \ 840520107fSSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \ 850520107fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; \ 860520107fSSatish Balay a->singlemalloc = 1; \ 870520107fSSatish Balay \ 880520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 8930770e4dSSatish Balay rmax = aimax[row] = aimax[row] + CHUNKSIZE; \ 900520107fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 910520107fSSatish Balay a->maxnz += CHUNKSIZE; \ 920520107fSSatish Balay a->reallocs++; \ 930520107fSSatish Balay } \ 940520107fSSatish Balay N = nrow++ - 1; a->nz++; \ 950520107fSSatish Balay /* shift up all the later entries in this row */ \ 960520107fSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 970520107fSSatish Balay rp[ii+1] = rp[ii]; \ 980520107fSSatish Balay ap[ii+1] = ap[ii]; \ 990520107fSSatish Balay } \ 100f5e9677aSSatish Balay rp[_i] = col1; \ 1010520107fSSatish Balay ap[_i] = value; \ 10230770e4dSSatish Balay a_noinsert: ; \ 1030520107fSSatish Balay ailen[row] = nrow; \ 1040520107fSSatish Balay } 1050a198c4cSBarry Smith 10630770e4dSSatish Balay #define MatSetValues_SeqAIJ_B_Private(row,col,value,addv) \ 10730770e4dSSatish Balay { \ 10830770e4dSSatish Balay \ 10930770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 11030770e4dSSatish Balay rmax = bimax[row]; nrow = bilen[row]; \ 11130770e4dSSatish Balay col1 = col - shift; \ 11230770e4dSSatish Balay \ 113ba4e3ef2SSatish Balay low = 0; high = nrow; \ 114ba4e3ef2SSatish Balay while (high-low > 5) { \ 115ba4e3ef2SSatish Balay t = (low+high)/2; \ 116ba4e3ef2SSatish Balay if (rp[t] > col) high = t; \ 117ba4e3ef2SSatish Balay else low = t; \ 118ba4e3ef2SSatish Balay } \ 11930770e4dSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 12030770e4dSSatish Balay if (rp[_i] > col1) break; \ 12130770e4dSSatish Balay if (rp[_i] == col1) { \ 12230770e4dSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 12330770e4dSSatish Balay else ap[_i] = value; \ 12430770e4dSSatish Balay goto b_noinsert; \ 12530770e4dSSatish Balay } \ 12630770e4dSSatish Balay } \ 12789280ab3SLois Curfman McInnes if (nonew == 1) goto b_noinsert; \ 12811ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 12930770e4dSSatish Balay if (nrow >= rmax) { \ 13030770e4dSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 13174c639caSSatish Balay int new_nz = bi[b->m] + CHUNKSIZE,len,*new_i,*new_j; \ 13230770e4dSSatish Balay Scalar *new_a; \ 13330770e4dSSatish Balay \ 13411ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 13589280ab3SLois Curfman McInnes \ 13630770e4dSSatish Balay /* malloc new storage space */ \ 13774c639caSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(b->m+1)*sizeof(int); \ 13830770e4dSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 13930770e4dSSatish Balay new_j = (int *) (new_a + new_nz); \ 14030770e4dSSatish Balay new_i = new_j + new_nz; \ 14130770e4dSSatish Balay \ 14230770e4dSSatish Balay /* copy over old data into new slots */ \ 14330770e4dSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = bi[ii];} \ 14474c639caSSatish Balay for ( ii=row+1; ii<b->m+1; ii++ ) {new_i[ii] = bi[ii]+CHUNKSIZE;} \ 14530770e4dSSatish Balay PetscMemcpy(new_j,bj,(bi[row]+nrow+shift)*sizeof(int)); \ 14630770e4dSSatish Balay len = (new_nz - CHUNKSIZE - bi[row] - nrow - shift); \ 14730770e4dSSatish Balay PetscMemcpy(new_j+bi[row]+shift+nrow+CHUNKSIZE,bj+bi[row]+shift+nrow, \ 14830770e4dSSatish Balay len*sizeof(int)); \ 14930770e4dSSatish Balay PetscMemcpy(new_a,ba,(bi[row]+nrow+shift)*sizeof(Scalar)); \ 15030770e4dSSatish Balay PetscMemcpy(new_a+bi[row]+shift+nrow+CHUNKSIZE,ba+bi[row]+shift+nrow, \ 15130770e4dSSatish Balay len*sizeof(Scalar)); \ 15230770e4dSSatish Balay /* free up old matrix storage */ \ 15330770e4dSSatish Balay \ 15474c639caSSatish Balay PetscFree(b->a); \ 15574c639caSSatish Balay if (!b->singlemalloc) {PetscFree(b->i);PetscFree(b->j);} \ 15674c639caSSatish Balay ba = b->a = new_a; bi = b->i = new_i; bj = b->j = new_j; \ 15774c639caSSatish Balay b->singlemalloc = 1; \ 15830770e4dSSatish Balay \ 15930770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 16030770e4dSSatish Balay rmax = bimax[row] = bimax[row] + CHUNKSIZE; \ 16174c639caSSatish Balay PLogObjectMemory(B,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 16274c639caSSatish Balay b->maxnz += CHUNKSIZE; \ 16374c639caSSatish Balay b->reallocs++; \ 16430770e4dSSatish Balay } \ 16574c639caSSatish Balay N = nrow++ - 1; b->nz++; \ 16630770e4dSSatish Balay /* shift up all the later entries in this row */ \ 16730770e4dSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 16830770e4dSSatish Balay rp[ii+1] = rp[ii]; \ 16930770e4dSSatish Balay ap[ii+1] = ap[ii]; \ 17030770e4dSSatish Balay } \ 17130770e4dSSatish Balay rp[_i] = col1; \ 17230770e4dSSatish Balay ap[_i] = value; \ 17330770e4dSSatish Balay b_noinsert: ; \ 17430770e4dSSatish Balay bilen[row] = nrow; \ 17530770e4dSSatish Balay } 17630770e4dSSatish Balay 1770520107fSSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode); 1785615d1e5SSatish Balay #undef __FUNC__ 1795615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIAIJ" 1808f6be9afSLois Curfman McInnes int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 1818a729477SBarry Smith { 18244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1834b0e389bSBarry Smith Scalar value; 1841eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 1851eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 186905e6a2fSBarry Smith int roworiented = aij->roworiented; 1878a729477SBarry Smith 1880520107fSSatish Balay /* Some Variables required in the macro */ 1894ee7247eSSatish Balay Mat A = aij->A; 1904ee7247eSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 19130770e4dSSatish Balay int *aimax = a->imax, *ai = a->i, *ailen = a->ilen,*aj = a->j; 19230770e4dSSatish Balay Scalar *aa = a->a; 19330770e4dSSatish Balay 19430770e4dSSatish Balay Mat B = aij->B; 19530770e4dSSatish Balay Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data; 19630770e4dSSatish Balay int *bimax = b->imax, *bi = b->i, *bilen = b->ilen,*bj = b->j; 19730770e4dSSatish Balay Scalar *ba = b->a; 19830770e4dSSatish Balay 199ba4e3ef2SSatish Balay int *rp,ii,nrow,_i,rmax, N, col1,low,high,t; 20030770e4dSSatish Balay int nonew = a->nonew,shift = a->indexshift; 20130770e4dSSatish Balay Scalar *ap; 2024ee7247eSSatish Balay 2033a40ed3dSBarry Smith PetscFunctionBegin; 2048a729477SBarry Smith for ( i=0; i<m; i++ ) { 2053a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 206e3372554SBarry Smith if (im[i] < 0) SETERRQ(1,0,"Negative row"); 207e3372554SBarry Smith if (im[i] >= aij->M) SETERRQ(1,0,"Row too large"); 2080a198c4cSBarry Smith #endif 2094b0e389bSBarry Smith if (im[i] >= rstart && im[i] < rend) { 2104b0e389bSBarry Smith row = im[i] - rstart; 2111eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 2124b0e389bSBarry Smith if (in[j] >= cstart && in[j] < cend){ 2134b0e389bSBarry Smith col = in[j] - cstart; 2144b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 21530770e4dSSatish Balay MatSetValues_SeqAIJ_A_Private(row,col,value,addv); 2160520107fSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2171eb62cbbSBarry Smith } 2183a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 219e3372554SBarry Smith else if (in[j] < 0) {SETERRQ(1,0,"Negative column");} 220e3372554SBarry Smith else if (in[j] >= aij->N) {SETERRQ(1,0,"Col too large");} 2210a198c4cSBarry Smith #endif 2221eb62cbbSBarry Smith else { 223227d817aSBarry Smith if (mat->was_assembled) { 224905e6a2fSBarry Smith if (!aij->colmap) { 225905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 226905e6a2fSBarry Smith } 227905e6a2fSBarry Smith col = aij->colmap[in[j]] - 1; 228ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 2292493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 2304b0e389bSBarry Smith col = in[j]; 2319bf004c3SSatish Balay /* Reinitialize the variables required by MatSetValues_SeqAIJ_B_Private() */ 232f9508a3cSSatish Balay B = aij->B; 233f9508a3cSSatish Balay b = (Mat_SeqAIJ *) B->data; 234f9508a3cSSatish Balay bimax = b->imax; bi = b->i; bilen = b->ilen; bj = b->j; 235f9508a3cSSatish Balay ba = b->a; 236d6dfbf8fSBarry Smith } 2379e25ed09SBarry Smith } 2384b0e389bSBarry Smith else col = in[j]; 2394b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 24030770e4dSSatish Balay MatSetValues_SeqAIJ_B_Private(row,col,value,addv); 24130770e4dSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2421eb62cbbSBarry Smith } 2431eb62cbbSBarry Smith } 2441eb62cbbSBarry Smith } 2451eb62cbbSBarry Smith else { 24690f02eecSBarry Smith if (roworiented && !aij->donotstash) { 2474b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr); 2484b0e389bSBarry Smith } 2494b0e389bSBarry Smith else { 25090f02eecSBarry Smith if (!aij->donotstash) { 2514b0e389bSBarry Smith row = im[i]; 2524b0e389bSBarry Smith for ( j=0; j<n; j++ ) { 2534b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr); 2544b0e389bSBarry Smith } 2554b0e389bSBarry Smith } 2561eb62cbbSBarry Smith } 2578a729477SBarry Smith } 25890f02eecSBarry Smith } 2593a40ed3dSBarry Smith PetscFunctionReturn(0); 2608a729477SBarry Smith } 2618a729477SBarry Smith 2625615d1e5SSatish Balay #undef __FUNC__ 2635615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIAIJ" 2648f6be9afSLois Curfman McInnes int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 265b49de8d1SLois Curfman McInnes { 266b49de8d1SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 267b49de8d1SLois Curfman McInnes int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 268b49de8d1SLois Curfman McInnes int cstart = aij->cstart, cend = aij->cend,row,col; 269b49de8d1SLois Curfman McInnes 2703a40ed3dSBarry Smith PetscFunctionBegin; 271b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 272e3372554SBarry Smith if (idxm[i] < 0) SETERRQ(1,0,"Negative row"); 273e3372554SBarry Smith if (idxm[i] >= aij->M) SETERRQ(1,0,"Row too large"); 274b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 275b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 276b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 277e3372554SBarry Smith if (idxn[j] < 0) SETERRQ(1,0,"Negative column"); 278e3372554SBarry Smith if (idxn[j] >= aij->N) SETERRQ(1,0,"Col too large"); 279b49de8d1SLois Curfman McInnes if (idxn[j] >= cstart && idxn[j] < cend){ 280b49de8d1SLois Curfman McInnes col = idxn[j] - cstart; 281b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 282b49de8d1SLois Curfman McInnes } 283b49de8d1SLois Curfman McInnes else { 284905e6a2fSBarry Smith if (!aij->colmap) { 285905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 286905e6a2fSBarry Smith } 287905e6a2fSBarry Smith col = aij->colmap[idxn[j]] - 1; 288e60e1c95SSatish Balay if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0; 289d9d09a02SSatish Balay else { 290b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 291b49de8d1SLois Curfman McInnes } 292b49de8d1SLois Curfman McInnes } 293b49de8d1SLois Curfman McInnes } 294d9d09a02SSatish Balay } 295b49de8d1SLois Curfman McInnes else { 296e3372554SBarry Smith SETERRQ(1,0,"Only local values currently supported"); 297b49de8d1SLois Curfman McInnes } 298b49de8d1SLois Curfman McInnes } 2993a40ed3dSBarry Smith PetscFunctionReturn(0); 300b49de8d1SLois Curfman McInnes } 301b49de8d1SLois Curfman McInnes 3025615d1e5SSatish Balay #undef __FUNC__ 3035615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ" 3048f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode) 3058a729477SBarry Smith { 30644a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 307d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 30817699dbbSLois Curfman McInnes int size = aij->size, *owners = aij->rowners; 30917699dbbSLois Curfman McInnes int rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr; 3101eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3116abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 3121eb62cbbSBarry Smith InsertMode addv; 3131eb62cbbSBarry Smith Scalar *rvalues,*svalues; 3141eb62cbbSBarry Smith 3153a40ed3dSBarry Smith PetscFunctionBegin; 3161eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 317ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 318dbb450caSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 319e3372554SBarry Smith SETERRQ(1,0,"Some processors inserted others added"); 3201eb62cbbSBarry Smith } 32147794344SBarry Smith mat->insertmode = addv; /* in case this processor had no cache */ 3221eb62cbbSBarry Smith 3231eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3240452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 325cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 3260452661fSBarry Smith owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 3271eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3281eb62cbbSBarry Smith idx = aij->stash.idx[i]; 32917699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 3301eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3311eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 3328a729477SBarry Smith } 3338a729477SBarry Smith } 3348a729477SBarry Smith } 33517699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 3361eb62cbbSBarry Smith 3371eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3380452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 339ca161407SBarry Smith ierr = MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 34017699dbbSLois Curfman McInnes nreceives = work[rank]; 341ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 34217699dbbSLois Curfman McInnes nmax = work[rank]; 3430452661fSBarry Smith PetscFree(work); 3441eb62cbbSBarry Smith 3451eb62cbbSBarry Smith /* post receives: 3461eb62cbbSBarry Smith 1) each message will consist of ordered pairs 3471eb62cbbSBarry Smith (global index,value) we store the global index as a double 348d6dfbf8fSBarry Smith to simplify the message passing. 3491eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 3501eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 3511eb62cbbSBarry Smith this is a lot of wasted space. 3521eb62cbbSBarry Smith 3531eb62cbbSBarry Smith 3541eb62cbbSBarry Smith This could be done better. 3551eb62cbbSBarry Smith */ 356ca161407SBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues); 357ca161407SBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 3581eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 359ca161407SBarry Smith ierr = MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 360ca161407SBarry Smith comm,recv_waits+i);CHKERRQ(ierr); 3611eb62cbbSBarry Smith } 3621eb62cbbSBarry Smith 3631eb62cbbSBarry Smith /* do sends: 3641eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3651eb62cbbSBarry Smith the ith processor 3661eb62cbbSBarry Smith */ 3670452661fSBarry Smith svalues = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 368ca161407SBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 3690452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 3701eb62cbbSBarry Smith starts[0] = 0; 37117699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3721eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3731eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 3741eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 3751eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 3761eb62cbbSBarry Smith } 3770452661fSBarry Smith PetscFree(owner); 3781eb62cbbSBarry Smith starts[0] = 0; 37917699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3801eb62cbbSBarry Smith count = 0; 38117699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 3821eb62cbbSBarry Smith if (procs[i]) { 383ca161407SBarry Smith ierr = MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3841eb62cbbSBarry Smith } 3851eb62cbbSBarry Smith } 3860452661fSBarry Smith PetscFree(starts); PetscFree(nprocs); 3871eb62cbbSBarry Smith 3881eb62cbbSBarry Smith /* Free cache space */ 38955908cccSLois Curfman McInnes PLogInfo(mat,"MatAssemblyBegin_MPIAIJ:Number of off-processor values %d\n",aij->stash.n); 39078b31e54SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 3911eb62cbbSBarry Smith 3921eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 3931eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 3941eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 3951eb62cbbSBarry Smith aij->rmax = nmax; 3961eb62cbbSBarry Smith 3973a40ed3dSBarry Smith PetscFunctionReturn(0); 3981eb62cbbSBarry Smith } 39944a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat); 4001eb62cbbSBarry Smith 4015615d1e5SSatish Balay #undef __FUNC__ 4025615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ" 4038f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode) 4041eb62cbbSBarry Smith { 40544a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 4061eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 407416022c9SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr; 408905e6a2fSBarry Smith int row,col,other_disassembled; 4091eb62cbbSBarry Smith Scalar *values,val; 41047794344SBarry Smith InsertMode addv = mat->insertmode; 4111eb62cbbSBarry Smith 4123a40ed3dSBarry Smith PetscFunctionBegin; 4131eb62cbbSBarry Smith /* wait on receives */ 4141eb62cbbSBarry Smith while (count) { 415ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4161eb62cbbSBarry Smith /* unpack receives into our local space */ 417d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 418ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_SCALAR,&n);CHKERRQ(ierr); 4191eb62cbbSBarry Smith n = n/3; 4201eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 421227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - aij->rstart; 422227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 4231eb62cbbSBarry Smith val = values[3*i+2]; 4241eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 4251eb62cbbSBarry Smith col -= aij->cstart; 4266fd7127cSSatish Balay ierr = MatSetValues(aij->A,1,&row,1,&col,&val,addv); CHKERRQ(ierr); 4273a40ed3dSBarry Smith } else { 42855a1b374SBarry Smith if (mat->was_assembled || mat->assembled) { 429905e6a2fSBarry Smith if (!aij->colmap) { 430905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat); CHKERRQ(ierr); 431905e6a2fSBarry Smith } 432905e6a2fSBarry Smith col = aij->colmap[col] - 1; 433ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 4342493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 435227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 436d6dfbf8fSBarry Smith } 4379e25ed09SBarry Smith } 4386fd7127cSSatish Balay ierr = MatSetValues(aij->B,1,&row,1,&col,&val,addv); CHKERRQ(ierr); 4391eb62cbbSBarry Smith } 4401eb62cbbSBarry Smith } 4411eb62cbbSBarry Smith count--; 4421eb62cbbSBarry Smith } 4430452661fSBarry Smith PetscFree(aij->recv_waits); PetscFree(aij->rvalues); 4441eb62cbbSBarry Smith 4451eb62cbbSBarry Smith /* wait on sends */ 4461eb62cbbSBarry Smith if (aij->nsends) { 4470a198c4cSBarry Smith send_status = (MPI_Status *) PetscMalloc(aij->nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 448ca161407SBarry Smith ierr = MPI_Waitall(aij->nsends,aij->send_waits,send_status);CHKERRQ(ierr); 4490452661fSBarry Smith PetscFree(send_status); 4501eb62cbbSBarry Smith } 4510452661fSBarry Smith PetscFree(aij->send_waits); PetscFree(aij->svalues); 4521eb62cbbSBarry Smith 45378b31e54SBarry Smith ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr); 45478b31e54SBarry Smith ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr); 4551eb62cbbSBarry Smith 4562493cbb0SBarry Smith /* determine if any processor has disassembled, if so we must 4572493cbb0SBarry Smith also disassemble ourselfs, in order that we may reassemble. */ 458ca161407SBarry Smith ierr = MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);CHKERRQ(ierr); 459227d817aSBarry Smith if (mat->was_assembled && !other_disassembled) { 4602493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 4612493cbb0SBarry Smith } 4622493cbb0SBarry Smith 4636d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 46478b31e54SBarry Smith ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr); 4655e42470aSBarry Smith } 46678b31e54SBarry Smith ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr); 46778b31e54SBarry Smith ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr); 4685e42470aSBarry Smith 4697a0afa10SBarry Smith if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;} 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 4718a729477SBarry Smith } 4728a729477SBarry Smith 4735615d1e5SSatish Balay #undef __FUNC__ 4745615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ" 4758f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A) 4761eb62cbbSBarry Smith { 47744a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 478dbd7a890SLois Curfman McInnes int ierr; 4793a40ed3dSBarry Smith 4803a40ed3dSBarry Smith PetscFunctionBegin; 48178b31e54SBarry Smith ierr = MatZeroEntries(l->A); CHKERRQ(ierr); 48278b31e54SBarry Smith ierr = MatZeroEntries(l->B); CHKERRQ(ierr); 4833a40ed3dSBarry Smith PetscFunctionReturn(0); 4841eb62cbbSBarry Smith } 4851eb62cbbSBarry Smith 4861eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 4871eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 4881eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 4891eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 4901eb62cbbSBarry Smith routine. 4911eb62cbbSBarry Smith */ 4925615d1e5SSatish Balay #undef __FUNC__ 4935615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ" 4948f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 4951eb62cbbSBarry Smith { 49644a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 49717699dbbSLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 4986abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 49917699dbbSLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 5005392566eSBarry Smith int *rvalues,tag = A->tag,count,base,slen,n,*source; 501d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 502d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 5031eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 5041eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 5051eb62cbbSBarry Smith IS istmp; 5061eb62cbbSBarry Smith 5073a40ed3dSBarry Smith PetscFunctionBegin; 50877c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 50978b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5101eb62cbbSBarry Smith 5111eb62cbbSBarry Smith /* first count number of contributors to each processor */ 5120452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 513cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 5140452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 5151eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5161eb62cbbSBarry Smith idx = rows[i]; 5171eb62cbbSBarry Smith found = 0; 51817699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 5191eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 5201eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 5211eb62cbbSBarry Smith } 5221eb62cbbSBarry Smith } 523e3372554SBarry Smith if (!found) SETERRQ(1,0,"Index out of range"); 5241eb62cbbSBarry Smith } 52517699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 5261eb62cbbSBarry Smith 5271eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 5280452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 529ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 53017699dbbSLois Curfman McInnes nrecvs = work[rank]; 531ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 53217699dbbSLois Curfman McInnes nmax = work[rank]; 5330452661fSBarry Smith PetscFree(work); 5341eb62cbbSBarry Smith 5351eb62cbbSBarry Smith /* post receives: */ 5363a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 537ca161407SBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 5381eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 539ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 5401eb62cbbSBarry Smith } 5411eb62cbbSBarry Smith 5421eb62cbbSBarry Smith /* do sends: 5431eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 5441eb62cbbSBarry Smith the ith processor 5451eb62cbbSBarry Smith */ 5460452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 5473a40ed3dSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 5480452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 5491eb62cbbSBarry Smith starts[0] = 0; 55017699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5511eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5521eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 5531eb62cbbSBarry Smith } 5541eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 5551eb62cbbSBarry Smith 5561eb62cbbSBarry Smith starts[0] = 0; 55717699dbbSLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5581eb62cbbSBarry Smith count = 0; 55917699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 5601eb62cbbSBarry Smith if (procs[i]) { 561ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 5621eb62cbbSBarry Smith } 5631eb62cbbSBarry Smith } 5640452661fSBarry Smith PetscFree(starts); 5651eb62cbbSBarry Smith 56617699dbbSLois Curfman McInnes base = owners[rank]; 5671eb62cbbSBarry Smith 5681eb62cbbSBarry Smith /* wait on receives */ 5690452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 5701eb62cbbSBarry Smith source = lens + nrecvs; 5711eb62cbbSBarry Smith count = nrecvs; slen = 0; 5721eb62cbbSBarry Smith while (count) { 573ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 5741eb62cbbSBarry Smith /* unpack receives into our local space */ 575ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 576d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 577d6dfbf8fSBarry Smith lens[imdex] = n; 5781eb62cbbSBarry Smith slen += n; 5791eb62cbbSBarry Smith count--; 5801eb62cbbSBarry Smith } 5810452661fSBarry Smith PetscFree(recv_waits); 5821eb62cbbSBarry Smith 5831eb62cbbSBarry Smith /* move the data into the send scatter */ 5840452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 5851eb62cbbSBarry Smith count = 0; 5861eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 5871eb62cbbSBarry Smith values = rvalues + i*nmax; 5881eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 5891eb62cbbSBarry Smith lrows[count++] = values[j] - base; 5901eb62cbbSBarry Smith } 5911eb62cbbSBarry Smith } 5920452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 5930452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 5941eb62cbbSBarry Smith 5951eb62cbbSBarry Smith /* actually zap the local rows */ 596029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 597464493b3SBarry Smith PLogObjectParent(A,istmp); 5980452661fSBarry Smith PetscFree(lrows); 59978b31e54SBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 60078b31e54SBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr); 60178b31e54SBarry Smith ierr = ISDestroy(istmp); CHKERRQ(ierr); 6021eb62cbbSBarry Smith 6031eb62cbbSBarry Smith /* wait on sends */ 6041eb62cbbSBarry Smith if (nsends) { 605ca161407SBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 606ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 6070452661fSBarry Smith PetscFree(send_status); 6081eb62cbbSBarry Smith } 6090452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 6101eb62cbbSBarry Smith 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6121eb62cbbSBarry Smith } 6131eb62cbbSBarry Smith 6145615d1e5SSatish Balay #undef __FUNC__ 6155615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ" 6168f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy) 6171eb62cbbSBarry Smith { 618416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 619fbd6ef76SBarry Smith int ierr,nt; 620416022c9SBarry Smith 6213a40ed3dSBarry Smith PetscFunctionBegin; 622a2ce50c7SBarry Smith ierr = VecGetLocalSize(xx,&nt); CHKERRQ(ierr); 623fbd6ef76SBarry Smith if (nt != a->n) { 624f40265daSLois Curfman McInnes SETERRQ(1,0,"Incompatible partition of A and xx"); 625fbd6ef76SBarry Smith } 62643a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 62738597bd4SSatish Balay ierr = (*a->A->ops.mult)(a->A,xx,yy); CHKERRQ(ierr); 62843a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 62938597bd4SSatish Balay ierr = (*a->B->ops.multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr); 6303a40ed3dSBarry Smith PetscFunctionReturn(0); 6311eb62cbbSBarry Smith } 6321eb62cbbSBarry Smith 6335615d1e5SSatish Balay #undef __FUNC__ 6345615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ" 6358f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 636da3a660dSBarry Smith { 637416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 638da3a660dSBarry Smith int ierr; 6393a40ed3dSBarry Smith 6403a40ed3dSBarry Smith PetscFunctionBegin; 64143a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 64238597bd4SSatish Balay ierr = (*a->A->ops.multadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 64343a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 64438597bd4SSatish Balay ierr = (*a->B->ops.multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr); 6453a40ed3dSBarry Smith PetscFunctionReturn(0); 646da3a660dSBarry Smith } 647da3a660dSBarry Smith 6485615d1e5SSatish Balay #undef __FUNC__ 6495615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ" 6508f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy) 651da3a660dSBarry Smith { 652416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 653da3a660dSBarry Smith int ierr; 654da3a660dSBarry Smith 6553a40ed3dSBarry Smith PetscFunctionBegin; 656da3a660dSBarry Smith /* do nondiagonal part */ 65738597bd4SSatish Balay ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 658da3a660dSBarry Smith /* send it on its way */ 659537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 660da3a660dSBarry Smith /* do local part */ 66138597bd4SSatish Balay ierr = (*a->A->ops.multtrans)(a->A,xx,yy); CHKERRQ(ierr); 662da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 663da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 664da3a660dSBarry Smith /* but is not perhaps always true. */ 665537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 6663a40ed3dSBarry Smith PetscFunctionReturn(0); 667da3a660dSBarry Smith } 668da3a660dSBarry Smith 6695615d1e5SSatish Balay #undef __FUNC__ 6705615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ" 6718f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 672da3a660dSBarry Smith { 673416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 674da3a660dSBarry Smith int ierr; 675da3a660dSBarry Smith 6763a40ed3dSBarry Smith PetscFunctionBegin; 677da3a660dSBarry Smith /* do nondiagonal part */ 67838597bd4SSatish Balay ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 679da3a660dSBarry Smith /* send it on its way */ 680537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 681da3a660dSBarry Smith /* do local part */ 68238597bd4SSatish Balay ierr = (*a->A->ops.multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 683da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 684da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 685da3a660dSBarry Smith /* but is not perhaps always true. */ 6860a198c4cSBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 6873a40ed3dSBarry Smith PetscFunctionReturn(0); 688da3a660dSBarry Smith } 689da3a660dSBarry Smith 6901eb62cbbSBarry Smith /* 6911eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 6921eb62cbbSBarry Smith diagonal block 6931eb62cbbSBarry Smith */ 6945615d1e5SSatish Balay #undef __FUNC__ 6955615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ" 6968f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v) 6971eb62cbbSBarry Smith { 6983a40ed3dSBarry Smith int ierr; 699416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 7003a40ed3dSBarry Smith 7013a40ed3dSBarry Smith PetscFunctionBegin; 7023a40ed3dSBarry Smith if (a->M != a->N) SETERRQ(1,0,"Supports only square matrix where A->A is diag block"); 7035baf8537SBarry Smith if (a->rstart != a->cstart || a->rend != a->cend) { 7043a40ed3dSBarry Smith SETERRQ(1,0,"row partition must equal col partition"); 7053a40ed3dSBarry Smith } 7063a40ed3dSBarry Smith ierr = MatGetDiagonal(a->A,v);CHKERRQ(ierr); 7073a40ed3dSBarry Smith PetscFunctionReturn(0); 7081eb62cbbSBarry Smith } 7091eb62cbbSBarry Smith 7105615d1e5SSatish Balay #undef __FUNC__ 7115615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ" 7128f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A) 713052efed2SBarry Smith { 714052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 715052efed2SBarry Smith int ierr; 7163a40ed3dSBarry Smith 7173a40ed3dSBarry Smith PetscFunctionBegin; 718052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 719052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 7203a40ed3dSBarry Smith PetscFunctionReturn(0); 721052efed2SBarry Smith } 722052efed2SBarry Smith 7235615d1e5SSatish Balay #undef __FUNC__ 724d4bb536fSBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" 7258f6be9afSLois Curfman McInnes int MatDestroy_MPIAIJ(PetscObject obj) 7261eb62cbbSBarry Smith { 7271eb62cbbSBarry Smith Mat mat = (Mat) obj; 72844a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 7291eb62cbbSBarry Smith int ierr; 73083e2fdc7SBarry Smith 7313a40ed3dSBarry Smith PetscFunctionBegin; 7323a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 7336e35fa57SLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N); 734a5a9c739SBarry Smith #endif 73583e2fdc7SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 7360452661fSBarry Smith PetscFree(aij->rowners); 73778b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 73878b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 7390452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 7400452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 7411eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 742a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 7437a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 7440452661fSBarry Smith PetscFree(aij); 745a5a9c739SBarry Smith PLogObjectDestroy(mat); 7460452661fSBarry Smith PetscHeaderDestroy(mat); 7473a40ed3dSBarry Smith PetscFunctionReturn(0); 7481eb62cbbSBarry Smith } 749ee50ffe9SBarry Smith 7505615d1e5SSatish Balay #undef __FUNC__ 751d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" 7528f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 7531eb62cbbSBarry Smith { 754416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 755416022c9SBarry Smith int ierr; 756416022c9SBarry Smith 7573a40ed3dSBarry Smith PetscFunctionBegin; 75817699dbbSLois Curfman McInnes if (aij->size == 1) { 759416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 760416022c9SBarry Smith } 761e3372554SBarry Smith else SETERRQ(1,0,"Only uniprocessor output supported"); 7623a40ed3dSBarry Smith PetscFunctionReturn(0); 763416022c9SBarry Smith } 764416022c9SBarry Smith 7655615d1e5SSatish Balay #undef __FUNC__ 766d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab" 7678f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer) 768416022c9SBarry Smith { 76944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 770dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 771a56f8943SBarry Smith int ierr, format,shift = C->indexshift,rank; 772d636dbe3SBarry Smith FILE *fd; 77319bcc07fSBarry Smith ViewerType vtype; 774416022c9SBarry Smith 7753a40ed3dSBarry Smith PetscFunctionBegin; 77619bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 77719bcc07fSBarry Smith if (vtype == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) { 77890ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 7790a198c4cSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 7804e220ebcSLois Curfman McInnes MatInfo info; 7814e220ebcSLois Curfman McInnes int flg; 782a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 78390ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 7844e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 78595e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 78677c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 78795e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 7884e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 78995e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 7904e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 7914e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->A,MAT_LOCAL,&info); 7924e220ebcSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used); 7934e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->B,MAT_LOCAL,&info); 7944e220ebcSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used); 795a56f8943SBarry Smith fflush(fd); 79677c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 797a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 7983a40ed3dSBarry Smith PetscFunctionReturn(0); 7993a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 8003a40ed3dSBarry Smith PetscFunctionReturn(0); 80108480c60SBarry Smith } 802416022c9SBarry Smith } 803416022c9SBarry Smith 80419bcc07fSBarry Smith if (vtype == DRAW_VIEWER) { 80519bcc07fSBarry Smith Draw draw; 80619bcc07fSBarry Smith PetscTruth isnull; 80719bcc07fSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 8083a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 80919bcc07fSBarry Smith } 81019bcc07fSBarry Smith 81119bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 81290ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 81377c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 814d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 81517699dbbSLois Curfman McInnes aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 8161eb62cbbSBarry Smith aij->cend); 81778b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 81878b31e54SBarry Smith ierr = MatView(aij->B,viewer); CHKERRQ(ierr); 819d13ab20cSBarry Smith fflush(fd); 82077c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 8213a40ed3dSBarry Smith } else { 822a56f8943SBarry Smith int size = aij->size; 823a56f8943SBarry Smith rank = aij->rank; 82417699dbbSLois Curfman McInnes if (size == 1) { 82578b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 8263a40ed3dSBarry Smith } else { 82795373324SBarry Smith /* assemble the entire matrix onto first processor. */ 82895373324SBarry Smith Mat A; 829ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 8302eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 83195373324SBarry Smith Scalar *a; 8322ee70a88SLois Curfman McInnes 83317699dbbSLois Curfman McInnes if (!rank) { 83455843e3eSBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr); 8353a40ed3dSBarry Smith } else { 83655843e3eSBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr); 83795373324SBarry Smith } 838464493b3SBarry Smith PLogObjectParent(mat,A); 839416022c9SBarry Smith 84095373324SBarry Smith /* copy over the A part */ 841ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 8422ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 84395373324SBarry Smith row = aij->rstart; 844dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 84595373324SBarry Smith for ( i=0; i<m; i++ ) { 846416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 84795373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 84895373324SBarry Smith } 8492ee70a88SLois Curfman McInnes aj = Aloc->j; 850dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 85195373324SBarry Smith 85295373324SBarry Smith /* copy over the B part */ 853ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 8542ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 85595373324SBarry Smith row = aij->rstart; 8560452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 857dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 85895373324SBarry Smith for ( i=0; i<m; i++ ) { 859416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 86095373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 86195373324SBarry Smith } 8620452661fSBarry Smith PetscFree(ct); 8636d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8646d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 86555843e3eSBarry Smith /* 86655843e3eSBarry Smith Everyone has to call to draw the matrix since the graphics waits are 86755843e3eSBarry Smith synchronized across all processors that share the Draw object 86855843e3eSBarry Smith */ 86955843e3eSBarry Smith if (!rank || vtype == DRAW_VIEWER) { 87078b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 87195373324SBarry Smith } 87278b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 87395373324SBarry Smith } 87495373324SBarry Smith } 8753a40ed3dSBarry Smith PetscFunctionReturn(0); 8761eb62cbbSBarry Smith } 8771eb62cbbSBarry Smith 8785615d1e5SSatish Balay #undef __FUNC__ 879d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ" 8808f6be9afSLois Curfman McInnes int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 881416022c9SBarry Smith { 882416022c9SBarry Smith Mat mat = (Mat) obj; 883416022c9SBarry Smith int ierr; 88419bcc07fSBarry Smith ViewerType vtype; 885416022c9SBarry Smith 8863a40ed3dSBarry Smith PetscFunctionBegin; 88719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 88819bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER || 88919bcc07fSBarry Smith vtype == DRAW_VIEWER || vtype == MATLAB_VIEWER) { 890d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 891416022c9SBarry Smith } 89219bcc07fSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 8933a40ed3dSBarry Smith ierr = MatView_MPIAIJ_Binary(mat,viewer);CHKERRQ(ierr); 894416022c9SBarry Smith } 8953a40ed3dSBarry Smith PetscFunctionReturn(0); 896416022c9SBarry Smith } 897416022c9SBarry Smith 8981eb62cbbSBarry Smith /* 8991eb62cbbSBarry Smith This has to provide several versions. 9001eb62cbbSBarry Smith 9011eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 9021eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 903d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 9041eb62cbbSBarry Smith */ 9055615d1e5SSatish Balay #undef __FUNC__ 9065615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ" 9078f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 908dbb450caSBarry Smith double fshift,int its,Vec xx) 9098a729477SBarry Smith { 91044a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 911d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 912ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 913c16cb8f2SBarry Smith Scalar *b,*x,*xs,*ls,d,*v,sum; 9146abc6512SBarry Smith int ierr,*idx, *diag; 915416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 9168a729477SBarry Smith 9173a40ed3dSBarry Smith PetscFunctionBegin; 918d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 919dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 920dbb450caSBarry Smith ls = ls + shift; 92183e2fdc7SBarry Smith if (!A->diag) {ierr = MatMarkDiag_SeqAIJ(AA); CHKERRQ(ierr);} 922d6dfbf8fSBarry Smith diag = A->diag; 923c16cb8f2SBarry Smith if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 924da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 9253a40ed3dSBarry Smith ierr = (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9263a40ed3dSBarry Smith PetscFunctionReturn(0); 927da3a660dSBarry Smith } 9283a40ed3dSBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 9293a40ed3dSBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 930d6dfbf8fSBarry Smith while (its--) { 931d6dfbf8fSBarry Smith /* go down through the rows */ 932d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 933d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 934dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 935dbb450caSBarry Smith v = A->a + A->i[i] + shift; 936d6dfbf8fSBarry Smith sum = b[i]; 937d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 938dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 939d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 940dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 941dbb450caSBarry Smith v = B->a + B->i[i] + shift; 942d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 94355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 944d6dfbf8fSBarry Smith } 945d6dfbf8fSBarry Smith /* come up through the rows */ 946d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 947d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 948dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 949dbb450caSBarry Smith v = A->a + A->i[i] + shift; 950d6dfbf8fSBarry Smith sum = b[i]; 951d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 952dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 953d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 954dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 955dbb450caSBarry Smith v = B->a + B->i[i] + shift; 956d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 95755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 958d6dfbf8fSBarry Smith } 959d6dfbf8fSBarry Smith } 9603a40ed3dSBarry Smith } else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 961da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 9623a40ed3dSBarry Smith ierr = (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9633a40ed3dSBarry Smith PetscFunctionReturn(0); 964da3a660dSBarry Smith } 9653a40ed3dSBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 9663a40ed3dSBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 967d6dfbf8fSBarry Smith while (its--) { 968d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 969d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 970dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 971dbb450caSBarry Smith v = A->a + A->i[i] + shift; 972d6dfbf8fSBarry Smith sum = b[i]; 973d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 974dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 975d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 976dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 977dbb450caSBarry Smith v = B->a + B->i[i] + shift; 978d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 97955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 980d6dfbf8fSBarry Smith } 981d6dfbf8fSBarry Smith } 9823a40ed3dSBarry Smith } else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 983da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 9843a40ed3dSBarry Smith ierr = (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9853a40ed3dSBarry Smith PetscFunctionReturn(0); 986da3a660dSBarry Smith } 98743a90d84SBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 98878b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 98943a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 99078b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 991d6dfbf8fSBarry Smith while (its--) { 992d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 993d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 994dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 995dbb450caSBarry Smith v = A->a + A->i[i] + shift; 996d6dfbf8fSBarry Smith sum = b[i]; 997d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 998dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 999d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 1000dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 1001dbb450caSBarry Smith v = B->a + B->i[i] + shift; 1002d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 100355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 1004d6dfbf8fSBarry Smith } 1005d6dfbf8fSBarry Smith } 10063a40ed3dSBarry Smith } else { 1007*b57c8cebSBarry Smith SETERRQ(1,0,"Parallel SOR not supported"); 1008c16cb8f2SBarry Smith } 10093a40ed3dSBarry Smith PetscFunctionReturn(0); 10108a729477SBarry Smith } 1011a66be287SLois Curfman McInnes 10125615d1e5SSatish Balay #undef __FUNC__ 1013d4bb536fSBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" 10148f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info) 1015a66be287SLois Curfman McInnes { 1016a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1017a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 10184e220ebcSLois Curfman McInnes int ierr; 10194e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 1020a66be287SLois Curfman McInnes 10213a40ed3dSBarry Smith PetscFunctionBegin; 10224e220ebcSLois Curfman McInnes info->block_size = 1.0; 10234e220ebcSLois Curfman McInnes ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr); 10244e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 10254e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 10264e220ebcSLois Curfman McInnes ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr); 10274e220ebcSLois Curfman McInnes isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded; 10284e220ebcSLois Curfman McInnes isend[3] += info->memory; isend[4] += info->mallocs; 1029a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 10304e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 10314e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 10324e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 10334e220ebcSLois Curfman McInnes info->memory = isend[3]; 10344e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 1035a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 1036ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm);CHKERRQ(ierr); 10374e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10384e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10394e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10404e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10414e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1042a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 1043ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm);CHKERRQ(ierr); 10444e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10454e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10464e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10474e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10484e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1049a66be287SLois Curfman McInnes } 10504e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 10514e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10524e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10538d700155SBarry Smith info->rows_global = (double)mat->M; 10548d700155SBarry Smith info->columns_global = (double)mat->N; 10558d700155SBarry Smith info->rows_local = (double)mat->m; 10568d700155SBarry Smith info->columns_local = (double)mat->N; 10574e220ebcSLois Curfman McInnes 10583a40ed3dSBarry Smith PetscFunctionReturn(0); 1059a66be287SLois Curfman McInnes } 1060a66be287SLois Curfman McInnes 10615615d1e5SSatish Balay #undef __FUNC__ 1062d4bb536fSBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" 10638f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op) 1064c74985f6SBarry Smith { 1065c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1066c74985f6SBarry Smith 10673a40ed3dSBarry Smith PetscFunctionBegin; 10686d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 10696d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 10706da5968aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 1071c2653b3dSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 107296854ed6SLois Curfman McInnes op == MAT_NEW_NONZERO_ALLOCATION_ERROR || 1073c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR) { 1074b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 1075b1fbbac0SLois Curfman McInnes MatSetOption(a->B,op); 1076b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 1077aeafbbfcSLois Curfman McInnes a->roworiented = 1; 1078c0bbcb79SLois Curfman McInnes MatSetOption(a->A,op); 1079c0bbcb79SLois Curfman McInnes MatSetOption(a->B,op); 1080b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 10816da5968aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 10826d4a8577SBarry Smith op == MAT_SYMMETRIC || 10836d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 10846d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 108594a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_MPIAIJ:Option ignored\n"); 10866d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) { 10874b0e389bSBarry Smith a->roworiented = 0; 10884b0e389bSBarry Smith MatSetOption(a->A,op); 10894b0e389bSBarry Smith MatSetOption(a->B,op); 10902b362799SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 109190f02eecSBarry Smith a->donotstash = 1; 10923a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS){ 10933a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 10943a40ed3dSBarry Smith } else { 10953a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 10963a40ed3dSBarry Smith } 10973a40ed3dSBarry Smith PetscFunctionReturn(0); 1098c74985f6SBarry Smith } 1099c74985f6SBarry Smith 11005615d1e5SSatish Balay #undef __FUNC__ 1101d4bb536fSBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" 11028f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1103c74985f6SBarry Smith { 110444a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 11053a40ed3dSBarry Smith 11063a40ed3dSBarry Smith PetscFunctionBegin; 11070752156aSBarry Smith if (m) *m = mat->M; 11080752156aSBarry Smith if (n) *n = mat->N; 11093a40ed3dSBarry Smith PetscFunctionReturn(0); 1110c74985f6SBarry Smith } 1111c74985f6SBarry Smith 11125615d1e5SSatish Balay #undef __FUNC__ 1113d4bb536fSBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" 11148f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1115c74985f6SBarry Smith { 111644a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 11173a40ed3dSBarry Smith 11183a40ed3dSBarry Smith PetscFunctionBegin; 11190752156aSBarry Smith if (m) *m = mat->m; 11200752156aSBarry Smith if (n) *n = mat->N; 11213a40ed3dSBarry Smith PetscFunctionReturn(0); 1122c74985f6SBarry Smith } 1123c74985f6SBarry Smith 11245615d1e5SSatish Balay #undef __FUNC__ 1125d4bb536fSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" 11268f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1127c74985f6SBarry Smith { 112844a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 11293a40ed3dSBarry Smith 11303a40ed3dSBarry Smith PetscFunctionBegin; 1131c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 11323a40ed3dSBarry Smith PetscFunctionReturn(0); 1133c74985f6SBarry Smith } 1134c74985f6SBarry Smith 11356d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11366d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11376d84be18SBarry Smith 11385615d1e5SSatish Balay #undef __FUNC__ 11395615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ" 11406d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 114139e00950SLois Curfman McInnes { 1142154123eaSLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 114370f0671dSBarry Smith Scalar *vworkA, *vworkB, **pvA, **pvB,*v_p; 1144154123eaSLois Curfman McInnes int i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart; 1145154123eaSLois Curfman McInnes int nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend; 114670f0671dSBarry Smith int *cmap, *idx_p; 114739e00950SLois Curfman McInnes 11483a40ed3dSBarry Smith PetscFunctionBegin; 1149e3372554SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,0,"Already active"); 11507a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11517a0afa10SBarry Smith 115270f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11537a0afa10SBarry Smith /* 11547a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11557a0afa10SBarry Smith */ 11567a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 1157c16cb8f2SBarry Smith int max = 1,m = mat->m,tmp; 1158c16cb8f2SBarry Smith for ( i=0; i<m; i++ ) { 11597a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11607a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11617a0afa10SBarry Smith } 11627a0afa10SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar))); 11637a0afa10SBarry Smith CHKPTRQ(mat->rowvalues); 11647a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11657a0afa10SBarry Smith } 11667a0afa10SBarry Smith 1167e3372554SBarry Smith if (row < rstart || row >= rend) SETERRQ(1,0,"Only local rows") 1168abc0e9e4SLois Curfman McInnes lrow = row - rstart; 116939e00950SLois Curfman McInnes 1170154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1171154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1172154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 117338597bd4SSatish Balay ierr = (*mat->A->ops.getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 117438597bd4SSatish Balay ierr = (*mat->B->ops.getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1175154123eaSLois Curfman McInnes nztot = nzA + nzB; 1176154123eaSLois Curfman McInnes 117770f0671dSBarry Smith cmap = mat->garray; 1178154123eaSLois Curfman McInnes if (v || idx) { 1179154123eaSLois Curfman McInnes if (nztot) { 1180154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 118170f0671dSBarry Smith int imark = -1; 1182154123eaSLois Curfman McInnes if (v) { 118370f0671dSBarry Smith *v = v_p = mat->rowvalues; 118439e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 118570f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1186154123eaSLois Curfman McInnes else break; 1187154123eaSLois Curfman McInnes } 1188154123eaSLois Curfman McInnes imark = i; 118970f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 119070f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1191154123eaSLois Curfman McInnes } 1192154123eaSLois Curfman McInnes if (idx) { 119370f0671dSBarry Smith *idx = idx_p = mat->rowindices; 119470f0671dSBarry Smith if (imark > -1) { 119570f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 119670f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 119770f0671dSBarry Smith } 119870f0671dSBarry Smith } else { 1199154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 120070f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1201154123eaSLois Curfman McInnes else break; 1202154123eaSLois Curfman McInnes } 1203154123eaSLois Curfman McInnes imark = i; 120470f0671dSBarry Smith } 120570f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 120670f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 120739e00950SLois Curfman McInnes } 120839e00950SLois Curfman McInnes } 12091ca473b0SSatish Balay else { 12101ca473b0SSatish Balay if (idx) *idx = 0; 12111ca473b0SSatish Balay if (v) *v = 0; 12121ca473b0SSatish Balay } 1213154123eaSLois Curfman McInnes } 121439e00950SLois Curfman McInnes *nz = nztot; 121538597bd4SSatish Balay ierr = (*mat->A->ops.restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 121638597bd4SSatish Balay ierr = (*mat->B->ops.restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 12173a40ed3dSBarry Smith PetscFunctionReturn(0); 121839e00950SLois Curfman McInnes } 121939e00950SLois Curfman McInnes 12205615d1e5SSatish Balay #undef __FUNC__ 1221d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" 12226d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 122339e00950SLois Curfman McInnes { 12247a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 12253a40ed3dSBarry Smith 12263a40ed3dSBarry Smith PetscFunctionBegin; 12277a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 1228e3372554SBarry Smith SETERRQ(1,0,"MatGetRow not called"); 12297a0afa10SBarry Smith } 12307a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 12313a40ed3dSBarry Smith PetscFunctionReturn(0); 123239e00950SLois Curfman McInnes } 123339e00950SLois Curfman McInnes 12345615d1e5SSatish Balay #undef __FUNC__ 12355615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ" 12368f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1237855ac2c5SLois Curfman McInnes { 1238855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1239ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1240416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1241416022c9SBarry Smith double sum = 0.0; 124204ca555eSLois Curfman McInnes Scalar *v; 124304ca555eSLois Curfman McInnes 12443a40ed3dSBarry Smith PetscFunctionBegin; 124517699dbbSLois Curfman McInnes if (aij->size == 1) { 124614183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 124737fa93a5SLois Curfman McInnes } else { 124804ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 124904ca555eSLois Curfman McInnes v = amat->a; 125004ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 12513a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 125204ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 125304ca555eSLois Curfman McInnes #else 125404ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 125504ca555eSLois Curfman McInnes #endif 125604ca555eSLois Curfman McInnes } 125704ca555eSLois Curfman McInnes v = bmat->a; 125804ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 12593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 126004ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 126104ca555eSLois Curfman McInnes #else 126204ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 126304ca555eSLois Curfman McInnes #endif 126404ca555eSLois Curfman McInnes } 1265ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr); 126604ca555eSLois Curfman McInnes *norm = sqrt(*norm); 12673a40ed3dSBarry Smith } else if (type == NORM_1) { /* max column norm */ 126804ca555eSLois Curfman McInnes double *tmp, *tmp2; 126904ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 1270758f045eSSatish Balay tmp = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp); 1271758f045eSSatish Balay tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2); 1272cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 127304ca555eSLois Curfman McInnes *norm = 0.0; 127404ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 127504ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1276579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 127704ca555eSLois Curfman McInnes } 127804ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 127904ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1280579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 128104ca555eSLois Curfman McInnes } 1282ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr); 128304ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 128404ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 128504ca555eSLois Curfman McInnes } 12860452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 12873a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 1288515d9167SLois Curfman McInnes double ntemp = 0.0; 128904ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1290dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 129104ca555eSLois Curfman McInnes sum = 0.0; 129204ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1293cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 129404ca555eSLois Curfman McInnes } 1295dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 129604ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1297cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 129804ca555eSLois Curfman McInnes } 1299515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 130004ca555eSLois Curfman McInnes } 1301ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);CHKERRQ(ierr); 1302ca161407SBarry Smith } else { 1303e3372554SBarry Smith SETERRQ(1,0,"No support for two norm"); 130404ca555eSLois Curfman McInnes } 130537fa93a5SLois Curfman McInnes } 13063a40ed3dSBarry Smith PetscFunctionReturn(0); 1307855ac2c5SLois Curfman McInnes } 1308855ac2c5SLois Curfman McInnes 13095615d1e5SSatish Balay #undef __FUNC__ 13105615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ" 13118f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1312b7c46309SBarry Smith { 1313b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1314dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1315416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1316b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 13173a40ed3dSBarry Smith Mat B; 1318b7c46309SBarry Smith Scalar *array; 1319b7c46309SBarry Smith 13203a40ed3dSBarry Smith PetscFunctionBegin; 1321d4bb536fSBarry Smith if (matout == PETSC_NULL && M != N) { 1322e3372554SBarry Smith SETERRQ(1,0,"Square matrix only for in-place"); 1323d4bb536fSBarry Smith } 1324d4bb536fSBarry Smith 1325d4bb536fSBarry Smith ierr = MatCreateMPIAIJ(A->comm,a->n,a->m,N,M,0,PETSC_NULL,0,PETSC_NULL,&B);CHKERRQ(ierr); 1326b7c46309SBarry Smith 1327b7c46309SBarry Smith /* copy over the A part */ 1328ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1329b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1330b7c46309SBarry Smith row = a->rstart; 1331dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1332b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1333416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1334b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1335b7c46309SBarry Smith } 1336b7c46309SBarry Smith aj = Aloc->j; 13374af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1338b7c46309SBarry Smith 1339b7c46309SBarry Smith /* copy over the B part */ 1340ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1341b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1342b7c46309SBarry Smith row = a->rstart; 13430452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1344dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1345b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1346416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1347b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1348b7c46309SBarry Smith } 13490452661fSBarry Smith PetscFree(ct); 13506d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13516d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13523638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13530de55854SLois Curfman McInnes *matout = B; 13540de55854SLois Curfman McInnes } else { 13550de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13560452661fSBarry Smith PetscFree(a->rowners); 13570de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13580de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 13590452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 13600452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13610de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1362a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13630452661fSBarry Smith PetscFree(a); 1364f09e8eb9SSatish Balay PetscMemcpy(A,B,sizeof(struct _p_Mat)); 13650452661fSBarry Smith PetscHeaderDestroy(B); 13660de55854SLois Curfman McInnes } 13673a40ed3dSBarry Smith PetscFunctionReturn(0); 1368b7c46309SBarry Smith } 1369b7c46309SBarry Smith 13705615d1e5SSatish Balay #undef __FUNC__ 13715615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ" 13724b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr) 1373a008b906SSatish Balay { 13744b967eb1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 13754b967eb1SSatish Balay Mat a = aij->A, b = aij->B; 1376a008b906SSatish Balay int ierr,s1,s2,s3; 1377a008b906SSatish Balay 13783a40ed3dSBarry Smith PetscFunctionBegin; 13794b967eb1SSatish Balay ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr); 13804b967eb1SSatish Balay if (rr) { 13814b967eb1SSatish Balay s3 = aij->n; 13824b967eb1SSatish Balay VecGetLocalSize_Fast(rr,s1); 1383e3372554SBarry Smith if (s1!=s3) SETERRQ(1,0,"right vector non-conforming local size"); 13844b967eb1SSatish Balay /* Overlap communication with computation. */ 138543a90d84SBarry Smith ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1386a008b906SSatish Balay } 13874b967eb1SSatish Balay if (ll) { 13884b967eb1SSatish Balay VecGetLocalSize_Fast(ll,s1); 1389e3372554SBarry Smith if (s1!=s2) SETERRQ(1,0,"left vector non-conforming local size"); 1390c351683dSSatish Balay ierr = (*b->ops.diagonalscale)(b,ll,0); CHKERRQ(ierr); 13914b967eb1SSatish Balay } 13924b967eb1SSatish Balay /* scale the diagonal block */ 1393c351683dSSatish Balay ierr = (*a->ops.diagonalscale)(a,ll,rr); CHKERRQ(ierr); 13944b967eb1SSatish Balay 13954b967eb1SSatish Balay if (rr) { 13964b967eb1SSatish Balay /* Do a scatter end and then right scale the off-diagonal block */ 139743a90d84SBarry Smith ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1398c351683dSSatish Balay ierr = (*b->ops.diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr); 13994b967eb1SSatish Balay } 14004b967eb1SSatish Balay 14013a40ed3dSBarry Smith PetscFunctionReturn(0); 1402a008b906SSatish Balay } 1403a008b906SSatish Balay 1404a008b906SSatish Balay 1405682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat); 14065615d1e5SSatish Balay #undef __FUNC__ 1407d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" 14088f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A) 1409682d7d0cSBarry Smith { 1410682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 14113a40ed3dSBarry Smith int ierr; 1412682d7d0cSBarry Smith 14133a40ed3dSBarry Smith PetscFunctionBegin; 14143a40ed3dSBarry Smith if (!a->rank) { 14153a40ed3dSBarry Smith ierr = MatPrintHelp_SeqAIJ(a->A);CHKERRQ(ierr); 14163a40ed3dSBarry Smith } 14173a40ed3dSBarry Smith PetscFunctionReturn(0); 1418682d7d0cSBarry Smith } 1419682d7d0cSBarry Smith 14205615d1e5SSatish Balay #undef __FUNC__ 1421d4bb536fSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" 14228f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs) 14235a838052SSatish Balay { 14243a40ed3dSBarry Smith PetscFunctionBegin; 14255a838052SSatish Balay *bs = 1; 14263a40ed3dSBarry Smith PetscFunctionReturn(0); 14275a838052SSatish Balay } 14285615d1e5SSatish Balay #undef __FUNC__ 1429d4bb536fSBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" 14308f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A) 1431bb5a7306SBarry Smith { 1432bb5a7306SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1433bb5a7306SBarry Smith int ierr; 14343a40ed3dSBarry Smith 14353a40ed3dSBarry Smith PetscFunctionBegin; 1436bb5a7306SBarry Smith ierr = MatSetUnfactored(a->A); CHKERRQ(ierr); 14373a40ed3dSBarry Smith PetscFunctionReturn(0); 1438bb5a7306SBarry Smith } 1439bb5a7306SBarry Smith 1440d4bb536fSBarry Smith #undef __FUNC__ 1441d4bb536fSBarry Smith #define __FUNC__ "MatEqual_MPIAIJ" 1442d4bb536fSBarry Smith int MatEqual_MPIAIJ(Mat A, Mat B, PetscTruth *flag) 1443d4bb536fSBarry Smith { 1444d4bb536fSBarry Smith Mat_MPIAIJ *matB = (Mat_MPIAIJ *) B->data,*matA = (Mat_MPIAIJ *) A->data; 1445d4bb536fSBarry Smith Mat a, b, c, d; 1446d4bb536fSBarry Smith PetscTruth flg; 1447d4bb536fSBarry Smith int ierr; 1448d4bb536fSBarry Smith 14493a40ed3dSBarry Smith PetscFunctionBegin; 1450d4bb536fSBarry Smith if (B->type != MATMPIAIJ) SETERRQ(1,0,"Matrices must be same type"); 1451d4bb536fSBarry Smith a = matA->A; b = matA->B; 1452d4bb536fSBarry Smith c = matB->A; d = matB->B; 1453d4bb536fSBarry Smith 1454d4bb536fSBarry Smith ierr = MatEqual(a, c, &flg); CHKERRQ(ierr); 1455d4bb536fSBarry Smith if (flg == PETSC_TRUE) { 1456d4bb536fSBarry Smith ierr = MatEqual(b, d, &flg); CHKERRQ(ierr); 1457d4bb536fSBarry Smith } 1458ca161407SBarry Smith ierr = MPI_Allreduce(&flg, flag, 1, MPI_INT, MPI_LAND, A->comm);CHKERRQ(ierr); 14593a40ed3dSBarry Smith PetscFunctionReturn(0); 1460d4bb536fSBarry Smith } 1461d4bb536fSBarry Smith 14628f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int); 14632f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 14640a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 14650a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **); 146600e6dbe6SBarry Smith extern int MatGetSubMatrix_MPIAIJ (Mat ,IS,IS,MatGetSubMatrixCall,Mat *); 146700e6dbe6SBarry Smith 14688a729477SBarry Smith /* -------------------------------------------------------------------*/ 14692ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ, 147039e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 147144a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 147244a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 147336ce4990SBarry Smith 0,0, 147436ce4990SBarry Smith 0,0, 147536ce4990SBarry Smith 0,0, 147644a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1477b7c46309SBarry Smith MatTranspose_MPIAIJ, 1478d4bb536fSBarry Smith MatGetInfo_MPIAIJ,MatEqual_MPIAIJ, 1479a008b906SSatish Balay MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ, 1480ee50ffe9SBarry Smith MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ, 14811eb62cbbSBarry Smith 0, 1482299609e3SLois Curfman McInnes MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ, 148336ce4990SBarry Smith 0,0,0,0, 1484d1710a03SLois Curfman McInnes MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ, 148536ce4990SBarry Smith 0,0, 148694a9d846SBarry Smith 0,0,MatConvertSameType_MPIAIJ,0,0, 1487b49de8d1SLois Curfman McInnes 0,0,0, 1488598137ffSSatish Balay MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0, 1489052efed2SBarry Smith MatPrintHelp_MPIAIJ, 14903b2fbd54SBarry Smith MatScale_MPIAIJ,0,0,0, 14910a198c4cSBarry Smith MatGetBlockSize_MPIAIJ,0,0,0,0, 149200e6dbe6SBarry Smith MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ, 1493ca161407SBarry Smith 0,0,MatGetSubMatrix_MPIAIJ}; 149436ce4990SBarry Smith 14958a729477SBarry Smith 14965615d1e5SSatish Balay #undef __FUNC__ 14975615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ" 14981987afe7SBarry Smith /*@C 1499ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 15003a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 15013a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 15023a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 15033a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 15048a729477SBarry Smith 15058a729477SBarry Smith Input Parameters: 15061eb62cbbSBarry Smith . comm - MPI communicator 15077d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 150892e8d321SLois Curfman McInnes This value should be the same as the local size used in creating the 150992e8d321SLois Curfman McInnes y vector for the matrix-vector product y = Ax. 15101a3896d6SBarry Smith . n - This value should be the same as the local size used in creating the 15111a3896d6SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 15121a3896d6SBarry Smith calculated if N is given) For square matrices n is almost always m. 15137d3e4905SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 151492e8d321SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1515ab693e5aSLois Curfman McInnes . d_nz - number of nonzeros per row in diagonal portion of local submatrix 1516ff756334SLois Curfman McInnes (same for all local rows) 15172bd5e0b2SLois Curfman McInnes . d_nzz - array containing the number of nonzeros in the various rows of the 151892e8d321SLois Curfman McInnes diagonal portion of the local submatrix (possibly different for each row) 15192bd5e0b2SLois Curfman McInnes or PETSC_NULL. You must leave room for the diagonal entry even if 15202bd5e0b2SLois Curfman McInnes it is zero. 15212bd5e0b2SLois Curfman McInnes . o_nz - number of nonzeros per row in the off-diagonal portion of local 1522ab693e5aSLois Curfman McInnes submatrix (same for all local rows). 15232bd5e0b2SLois Curfman McInnes . o_nzz - array containing the number of nonzeros in the various rows of the 15242bd5e0b2SLois Curfman McInnes off-diagonal portion of the local submatrix (possibly different for 15252bd5e0b2SLois Curfman McInnes each row) or PETSC_NULL. 15268a729477SBarry Smith 1527ff756334SLois Curfman McInnes Output Parameter: 152844cd7ae7SLois Curfman McInnes . A - the matrix 15298a729477SBarry Smith 1530ff756334SLois Curfman McInnes Notes: 1531ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1532ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 15330002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 15340002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1535ff756334SLois Curfman McInnes 1536ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1537ff756334SLois Curfman McInnes (possibly both). 1538ff756334SLois Curfman McInnes 15395511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 15405511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 15415511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 15425511cfe3SLois Curfman McInnes 15435511cfe3SLois Curfman McInnes Options Database Keys: 15445511cfe3SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 15456e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 15466e62573dSLois Curfman McInnes $ (max limit=5) 15476e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 15486e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 15496e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 15505511cfe3SLois Curfman McInnes 1551e0245417SLois Curfman McInnes Storage Information: 1552e0245417SLois Curfman McInnes For a square global matrix we define each processor's diagonal portion 1553e0245417SLois Curfman McInnes to be its local rows and the corresponding columns (a square submatrix); 1554e0245417SLois Curfman McInnes each processor's off-diagonal portion encompasses the remainder of the 1555e0245417SLois Curfman McInnes local matrix (a rectangular submatrix). 1556e0245417SLois Curfman McInnes 1557e0245417SLois Curfman McInnes The user can specify preallocated storage for the diagonal part of 15585ace5be8SLois Curfman McInnes the local submatrix with either d_nz or d_nnz (not both). Set 15595ace5be8SLois Curfman McInnes d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic 15605ace5be8SLois Curfman McInnes memory allocation. Likewise, specify preallocated storage for the 15615ace5be8SLois Curfman McInnes off-diagonal part of the local submatrix with o_nz or o_nnz (not both). 1562ff756334SLois Curfman McInnes 15635511cfe3SLois Curfman McInnes Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In 15645511cfe3SLois Curfman McInnes the figure below we depict these three local rows and all columns (0-11). 15652191d07cSBarry Smith 1566b810aeb4SBarry Smith $ 0 1 2 3 4 5 6 7 8 9 10 11 1567b810aeb4SBarry Smith $ ------------------- 15685511cfe3SLois Curfman McInnes $ row 3 | o o o d d d o o o o o o 15695511cfe3SLois Curfman McInnes $ row 4 | o o o d d d o o o o o o 15705511cfe3SLois Curfman McInnes $ row 5 | o o o d d d o o o o o o 15715511cfe3SLois Curfman McInnes $ ------------------- 1572b810aeb4SBarry Smith $ 1573b810aeb4SBarry Smith 15745511cfe3SLois Curfman McInnes Thus, any entries in the d locations are stored in the d (diagonal) 15755511cfe3SLois Curfman McInnes submatrix, and any entries in the o locations are stored in the 15765511cfe3SLois Curfman McInnes o (off-diagonal) submatrix. Note that the d and the o submatrices are 15775511cfe3SLois Curfman McInnes stored simply in the MATSEQAIJ format for compressed row storage. 15785511cfe3SLois Curfman McInnes 15795511cfe3SLois Curfman McInnes Now d_nz should indicate the number of nonzeros per row in the d matrix, 15805511cfe3SLois Curfman McInnes and o_nz should indicate the number of nonzeros per row in the o matrix. 15815511cfe3SLois Curfman McInnes In general, for PDE problems in which most nonzeros are near the diagonal, 15823d323bbdSBarry Smith one expects d_nz >> o_nz. For large problems you MUST preallocate memory 158392e8d321SLois Curfman McInnes or you will get TERRIBLE performance; see the users' manual chapter on 15846da5968aSLois Curfman McInnes matrices. 15853a511b96SLois Curfman McInnes 1586dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1587ff756334SLois Curfman McInnes 1588fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 15898a729477SBarry Smith @*/ 15901eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 159144cd7ae7SLois Curfman McInnes int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A) 15928a729477SBarry Smith { 159344cd7ae7SLois Curfman McInnes Mat B; 159444cd7ae7SLois Curfman McInnes Mat_MPIAIJ *b; 159536ce4990SBarry Smith int ierr, i,sum[2],work[2],size; 1596416022c9SBarry Smith 15973a40ed3dSBarry Smith PetscFunctionBegin; 15983914022bSBarry Smith MPI_Comm_size(comm,&size); 15993914022bSBarry Smith if (size == 1) { 16003914022bSBarry Smith if (M == PETSC_DECIDE) M = m; 16013914022bSBarry Smith if (N == PETSC_DECIDE) N = n; 16023914022bSBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,d_nz,d_nnz,A); CHKERRQ(ierr); 16033a40ed3dSBarry Smith PetscFunctionReturn(0); 16043914022bSBarry Smith } 16053914022bSBarry Smith 160644cd7ae7SLois Curfman McInnes *A = 0; 1607d4bb536fSBarry Smith PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATMPIAIJ,comm,MatDestroy,MatView); 160844cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 160944cd7ae7SLois Curfman McInnes B->data = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b); 161044cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_MPIAIJ)); 161144cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 161244cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_MPIAIJ; 161344cd7ae7SLois Curfman McInnes B->view = MatView_MPIAIJ; 161444cd7ae7SLois Curfman McInnes B->factor = 0; 161544cd7ae7SLois Curfman McInnes B->assembled = PETSC_FALSE; 161690f02eecSBarry Smith B->mapping = 0; 1617d6dfbf8fSBarry Smith 161847794344SBarry Smith B->insertmode = NOT_SET_VALUES; 16199eb4d147SSatish Balay b->size = size; 162044cd7ae7SLois Curfman McInnes MPI_Comm_rank(comm,&b->rank); 16211eb62cbbSBarry Smith 16223a40ed3dSBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) { 1623e3372554SBarry Smith SETERRQ(1,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz"); 16243a40ed3dSBarry Smith } 16251987afe7SBarry Smith 1626dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE || N == PETSC_DECIDE) { 16271eb62cbbSBarry Smith work[0] = m; work[1] = n; 1628ca161407SBarry Smith ierr = MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm );CHKERRQ(ierr); 1629dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE) M = sum[0]; 1630dbd7a890SLois Curfman McInnes if (N == PETSC_DECIDE) N = sum[1]; 16311eb62cbbSBarry Smith } 163244cd7ae7SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);} 163344cd7ae7SLois Curfman McInnes if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);} 163444cd7ae7SLois Curfman McInnes b->m = m; B->m = m; 163544cd7ae7SLois Curfman McInnes b->n = n; B->n = n; 163644cd7ae7SLois Curfman McInnes b->N = N; B->N = N; 163744cd7ae7SLois Curfman McInnes b->M = M; B->M = M; 16381eb62cbbSBarry Smith 16391eb62cbbSBarry Smith /* build local table of row and column ownerships */ 164044cd7ae7SLois Curfman McInnes b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners); 1641f09e8eb9SSatish Balay PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1642603f58a4SSatish Balay b->cowners = b->rowners + b->size + 2; 1643ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 164444cd7ae7SLois Curfman McInnes b->rowners[0] = 0; 164544cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 164644cd7ae7SLois Curfman McInnes b->rowners[i] += b->rowners[i-1]; 16478a729477SBarry Smith } 164844cd7ae7SLois Curfman McInnes b->rstart = b->rowners[b->rank]; 164944cd7ae7SLois Curfman McInnes b->rend = b->rowners[b->rank+1]; 1650ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 165144cd7ae7SLois Curfman McInnes b->cowners[0] = 0; 165244cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 165344cd7ae7SLois Curfman McInnes b->cowners[i] += b->cowners[i-1]; 16548a729477SBarry Smith } 165544cd7ae7SLois Curfman McInnes b->cstart = b->cowners[b->rank]; 165644cd7ae7SLois Curfman McInnes b->cend = b->cowners[b->rank+1]; 16578a729477SBarry Smith 16585ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 1659029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr); 166044cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->A); 16617b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 1662029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr); 166344cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->B); 16648a729477SBarry Smith 16651eb62cbbSBarry Smith /* build cache for off array entries formed */ 166644cd7ae7SLois Curfman McInnes ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr); 166790f02eecSBarry Smith b->donotstash = 0; 166844cd7ae7SLois Curfman McInnes b->colmap = 0; 166944cd7ae7SLois Curfman McInnes b->garray = 0; 167044cd7ae7SLois Curfman McInnes b->roworiented = 1; 16718a729477SBarry Smith 16721eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 167344cd7ae7SLois Curfman McInnes b->lvec = 0; 167444cd7ae7SLois Curfman McInnes b->Mvctx = 0; 16758a729477SBarry Smith 16767a0afa10SBarry Smith /* stuff for MatGetRow() */ 167744cd7ae7SLois Curfman McInnes b->rowindices = 0; 167844cd7ae7SLois Curfman McInnes b->rowvalues = 0; 167944cd7ae7SLois Curfman McInnes b->getrowactive = PETSC_FALSE; 16807a0afa10SBarry Smith 168144cd7ae7SLois Curfman McInnes *A = B; 16823a40ed3dSBarry Smith PetscFunctionReturn(0); 1683d6dfbf8fSBarry Smith } 1684c74985f6SBarry Smith 16855615d1e5SSatish Balay #undef __FUNC__ 16865615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ" 16878f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues) 1688d6dfbf8fSBarry Smith { 1689d6dfbf8fSBarry Smith Mat mat; 1690416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 1691a1b97e82SLois Curfman McInnes int ierr, len=0, flg; 1692d6dfbf8fSBarry Smith 16933a40ed3dSBarry Smith PetscFunctionBegin; 1694416022c9SBarry Smith *newmat = 0; 1695d4bb536fSBarry Smith PetscHeaderCreate(mat,_p_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm,MatDestroy,MatView); 1696a5a9c739SBarry Smith PLogObjectCreate(mat); 16970452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1698416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 169944a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 170044a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1701d6dfbf8fSBarry Smith mat->factor = matin->factor; 1702c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1703d6dfbf8fSBarry Smith 170444cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 170544cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 170644cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 170744cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 1708d6dfbf8fSBarry Smith 1709416022c9SBarry Smith a->rstart = oldmat->rstart; 1710416022c9SBarry Smith a->rend = oldmat->rend; 1711416022c9SBarry Smith a->cstart = oldmat->cstart; 1712416022c9SBarry Smith a->cend = oldmat->cend; 171317699dbbSLois Curfman McInnes a->size = oldmat->size; 171417699dbbSLois Curfman McInnes a->rank = oldmat->rank; 171547794344SBarry Smith mat->insertmode = NOT_SET_VALUES; 1716bcd2baecSBarry Smith a->rowvalues = 0; 1717bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1718d6dfbf8fSBarry Smith 1719603f58a4SSatish Balay a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1720f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1721603f58a4SSatish Balay a->cowners = a->rowners + a->size + 2; 1722603f58a4SSatish Balay PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int)); 1723416022c9SBarry Smith ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 17242ee70a88SLois Curfman McInnes if (oldmat->colmap) { 17250452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1726416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1727416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1728416022c9SBarry Smith } else a->colmap = 0; 17293f41c07dSBarry Smith if (oldmat->garray) { 17303f41c07dSBarry Smith len = ((Mat_SeqAIJ *) (oldmat->B->data))->n; 17313f41c07dSBarry Smith a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray); 1732464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 17333f41c07dSBarry Smith if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1734416022c9SBarry Smith } else a->garray = 0; 1735d6dfbf8fSBarry Smith 1736416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1737416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1738a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1739416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 1740416022c9SBarry Smith ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 1741416022c9SBarry Smith PLogObjectParent(mat,a->A); 1742416022c9SBarry Smith ierr = MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr); 1743416022c9SBarry Smith PLogObjectParent(mat,a->B); 17445dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 174525cdf11fSBarry Smith if (flg) { 1746682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1747682d7d0cSBarry Smith } 17488a729477SBarry Smith *newmat = mat; 17493a40ed3dSBarry Smith PetscFunctionReturn(0); 17508a729477SBarry Smith } 1751416022c9SBarry Smith 175277c4ece6SBarry Smith #include "sys.h" 1753416022c9SBarry Smith 17545615d1e5SSatish Balay #undef __FUNC__ 17555615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ" 175619bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat) 1757416022c9SBarry Smith { 1758d65a2f8fSBarry Smith Mat A; 1759d65a2f8fSBarry Smith Scalar *vals,*svals; 176019bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 1761416022c9SBarry Smith MPI_Status status; 17623a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 176317699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1764d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 176519bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 1766416022c9SBarry Smith 17673a40ed3dSBarry Smith PetscFunctionBegin; 176817699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 176917699dbbSLois Curfman McInnes if (!rank) { 177090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 17710752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr); 1772e3372554SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object"); 1773d64ed03dSBarry Smith if (header[3] < 0) { 1774d64ed03dSBarry Smith SETERRQ(1,1,"Matrix stored in special format on disk, cannot load as MPIAIJ"); 1775d64ed03dSBarry Smith } 17766c5fab8fSBarry Smith } 17776c5fab8fSBarry Smith 1778d64ed03dSBarry Smith 1779ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 1780416022c9SBarry Smith M = header[1]; N = header[2]; 1781416022c9SBarry Smith /* determine ownership of all rows */ 178217699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 17830452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1784ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1785416022c9SBarry Smith rowners[0] = 0; 178617699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1787416022c9SBarry Smith rowners[i] += rowners[i-1]; 1788416022c9SBarry Smith } 178917699dbbSLois Curfman McInnes rstart = rowners[rank]; 179017699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1791416022c9SBarry Smith 1792416022c9SBarry Smith /* distribute row lengths to all processors */ 17930452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1794416022c9SBarry Smith offlens = ourlens + (rend-rstart); 179517699dbbSLois Curfman McInnes if (!rank) { 17960452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 17970752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 17980452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 179917699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1800ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 18010452661fSBarry Smith PetscFree(sndcounts); 18023a40ed3dSBarry Smith } else { 1803ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 1804416022c9SBarry Smith } 1805416022c9SBarry Smith 180617699dbbSLois Curfman McInnes if (!rank) { 1807416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 18080452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1809cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 181017699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1811416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1812416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1813416022c9SBarry Smith } 1814416022c9SBarry Smith } 18150452661fSBarry Smith PetscFree(rowlengths); 1816416022c9SBarry Smith 1817416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1818416022c9SBarry Smith maxnz = 0; 181917699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 18200452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1821416022c9SBarry Smith } 18220452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1823416022c9SBarry Smith 1824416022c9SBarry Smith /* read in my part of the matrix column indices */ 1825416022c9SBarry Smith nz = procsnz[0]; 18260452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 18270752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr); 1828d65a2f8fSBarry Smith 1829d65a2f8fSBarry Smith /* read in every one elses and ship off */ 183017699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1831d65a2f8fSBarry Smith nz = procsnz[i]; 18320752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 1833ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 1834d65a2f8fSBarry Smith } 18350452661fSBarry Smith PetscFree(cols); 18363a40ed3dSBarry Smith } else { 1837416022c9SBarry Smith /* determine buffer space needed for message */ 1838416022c9SBarry Smith nz = 0; 1839416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1840416022c9SBarry Smith nz += ourlens[i]; 1841416022c9SBarry Smith } 18420452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1843416022c9SBarry Smith 1844416022c9SBarry Smith /* receive message of column indices*/ 1845ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1846ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1847e3372554SBarry Smith if (maxnz != nz) SETERRQ(1,0,"something is wrong with file"); 1848416022c9SBarry Smith } 1849416022c9SBarry Smith 1850416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 1851cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 1852416022c9SBarry Smith jj = 0; 1853416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1854416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 1855d65a2f8fSBarry Smith if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1856416022c9SBarry Smith jj++; 1857416022c9SBarry Smith } 1858416022c9SBarry Smith } 1859d65a2f8fSBarry Smith 1860d65a2f8fSBarry Smith /* create our matrix */ 1861416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1862416022c9SBarry Smith ourlens[i] -= offlens[i]; 1863416022c9SBarry Smith } 1864d65a2f8fSBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 1865d65a2f8fSBarry Smith A = *newmat; 18666d4a8577SBarry Smith MatSetOption(A,MAT_COLUMNS_SORTED); 1867d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1868d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 1869d65a2f8fSBarry Smith } 1870416022c9SBarry Smith 187117699dbbSLois Curfman McInnes if (!rank) { 18720452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1873416022c9SBarry Smith 1874416022c9SBarry Smith /* read in my part of the matrix numerical values */ 1875416022c9SBarry Smith nz = procsnz[0]; 18760752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1877d65a2f8fSBarry Smith 1878d65a2f8fSBarry Smith /* insert into matrix */ 1879d65a2f8fSBarry Smith jj = rstart; 1880d65a2f8fSBarry Smith smycols = mycols; 1881d65a2f8fSBarry Smith svals = vals; 1882d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1883d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1884d65a2f8fSBarry Smith smycols += ourlens[i]; 1885d65a2f8fSBarry Smith svals += ourlens[i]; 1886d65a2f8fSBarry Smith jj++; 1887416022c9SBarry Smith } 1888416022c9SBarry Smith 1889d65a2f8fSBarry Smith /* read in other processors and ship out */ 189017699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1891416022c9SBarry Smith nz = procsnz[i]; 18920752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1893ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 1894416022c9SBarry Smith } 18950452661fSBarry Smith PetscFree(procsnz); 18963a40ed3dSBarry Smith } else { 1897d65a2f8fSBarry Smith /* receive numeric values */ 18980452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1899416022c9SBarry Smith 1900d65a2f8fSBarry Smith /* receive message of values*/ 1901ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1902ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1903e3372554SBarry Smith if (maxnz != nz) SETERRQ(1,0,"something is wrong with file"); 1904d65a2f8fSBarry Smith 1905d65a2f8fSBarry Smith /* insert into matrix */ 1906d65a2f8fSBarry Smith jj = rstart; 1907d65a2f8fSBarry Smith smycols = mycols; 1908d65a2f8fSBarry Smith svals = vals; 1909d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1910d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1911d65a2f8fSBarry Smith smycols += ourlens[i]; 1912d65a2f8fSBarry Smith svals += ourlens[i]; 1913d65a2f8fSBarry Smith jj++; 1914d65a2f8fSBarry Smith } 1915d65a2f8fSBarry Smith } 19160452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1917d65a2f8fSBarry Smith 19186d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 19196d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 19203a40ed3dSBarry Smith PetscFunctionReturn(0); 1921416022c9SBarry Smith } 1922a0ff6018SBarry Smith 192329da9460SBarry Smith #undef __FUNC__ 192429da9460SBarry Smith #define __FUNC__ "MatGetSubMatrix_MPIAIJ" 1925a0ff6018SBarry Smith /* 192629da9460SBarry Smith Not great since it makes two copies of the submatrix, first an SeqAIJ 192729da9460SBarry Smith in local and then by concatenating the local matrices the end result. 192829da9460SBarry Smith Writing it directly would be much like MatGetSubMatrices_MPIAIJ() 1929a0ff6018SBarry Smith */ 1930a0ff6018SBarry Smith int MatGetSubMatrix_MPIAIJ(Mat mat,IS isrow,IS iscol,MatGetSubMatrixCall call,Mat *newmat) 1931a0ff6018SBarry Smith { 193200e6dbe6SBarry Smith int ierr, i, m,n,rstart,row,rend,nz,*cwork,size,rank,j; 1933a0ff6018SBarry Smith Mat *local,M; 193400e6dbe6SBarry Smith Scalar *vwork,*aa; 193500e6dbe6SBarry Smith MPI_Comm comm = mat->comm; 193600e6dbe6SBarry Smith Mat_SeqAIJ *aij; 193700e6dbe6SBarry Smith int *ii, *jj,nlocal,*dlens,*olens,dlen,olen,jend; 1938a0ff6018SBarry Smith 1939a0ff6018SBarry Smith PetscFunctionBegin; 194000e6dbe6SBarry Smith MPI_Comm_rank(comm,&rank); 194100e6dbe6SBarry Smith MPI_Comm_size(comm,&size); 194200e6dbe6SBarry Smith 1943a0ff6018SBarry Smith ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr); 1944a0ff6018SBarry Smith 1945a0ff6018SBarry Smith /* 1946a0ff6018SBarry Smith m - number of local rows 1947a0ff6018SBarry Smith n - number of columns (same on all processors) 1948a0ff6018SBarry Smith rstart - first row in new global matrix generated 1949a0ff6018SBarry Smith */ 1950a0ff6018SBarry Smith ierr = MatGetSize(*local,&m,&n);CHKERRQ(ierr); 1951a0ff6018SBarry Smith if (call == MAT_INITIAL_MATRIX) { 195200e6dbe6SBarry Smith aij = (Mat_SeqAIJ *) (*local)->data; 195300e6dbe6SBarry Smith if (aij->indexshift) SETERRQ(1,1,"No support for index shifted matrix"); 195400e6dbe6SBarry Smith ii = aij->i; 195500e6dbe6SBarry Smith jj = aij->j; 195600e6dbe6SBarry Smith 1957a0ff6018SBarry Smith /* 195800e6dbe6SBarry Smith Determine the number of non-zeros in the diagonal and off-diagonal 195900e6dbe6SBarry Smith portions of the matrix in order to do correct preallocation 1960a0ff6018SBarry Smith */ 196100e6dbe6SBarry Smith 196200e6dbe6SBarry Smith /* first get start and end of "diagonal" columns */ 196300e6dbe6SBarry Smith nlocal = n/size + ((n % size) > rank); 1964ca161407SBarry Smith ierr = MPI_Scan(&nlocal,&rend,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 196500e6dbe6SBarry Smith rstart = rend - nlocal; 196600e6dbe6SBarry Smith 196700e6dbe6SBarry Smith /* next, compute all the lengths */ 196800e6dbe6SBarry Smith dlens = (int *) PetscMalloc( (2*m+1)*sizeof(int) );CHKPTRQ(dlens); 196900e6dbe6SBarry Smith olens = dlens + m; 197000e6dbe6SBarry Smith for ( i=0; i<m; i++ ) { 197100e6dbe6SBarry Smith jend = ii[i+1] - ii[i]; 197200e6dbe6SBarry Smith olen = 0; 197300e6dbe6SBarry Smith dlen = 0; 197400e6dbe6SBarry Smith for ( j=0; j<jend; j++ ) { 197500e6dbe6SBarry Smith if ( *jj < rstart || *jj >= rend) olen++; 197600e6dbe6SBarry Smith else dlen++; 197700e6dbe6SBarry Smith jj++; 197800e6dbe6SBarry Smith } 197900e6dbe6SBarry Smith olens[i] = olen; 198000e6dbe6SBarry Smith dlens[i] = dlen; 198100e6dbe6SBarry Smith } 198200e6dbe6SBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,PETSC_DECIDE,n,0,dlens,0,olens,&M);CHKERRQ(ierr); 198300e6dbe6SBarry Smith PetscFree(dlens); 1984a0ff6018SBarry Smith } else { 1985a0ff6018SBarry Smith int ml,nl; 1986a0ff6018SBarry Smith 1987a0ff6018SBarry Smith M = *newmat; 1988a0ff6018SBarry Smith ierr = MatGetLocalSize(M,&ml,&nl);CHKERRQ(ierr); 1989a0ff6018SBarry Smith if (ml != m) SETERRQ(1,1,"Previous matrix must be same size/layout as request"); 1990a0ff6018SBarry Smith ierr = MatZeroEntries(M);CHKERRQ(ierr); 1991a0ff6018SBarry Smith } 1992a0ff6018SBarry Smith ierr = MatGetOwnershipRange(M,&rstart,&rend); CHKERRQ(ierr); 199300e6dbe6SBarry Smith aij = (Mat_SeqAIJ *) (*local)->data; 199400e6dbe6SBarry Smith if (aij->indexshift) SETERRQ(1,1,"No support for index shifted matrix"); 199500e6dbe6SBarry Smith ii = aij->i; 199600e6dbe6SBarry Smith jj = aij->j; 199700e6dbe6SBarry Smith aa = aij->a; 1998a0ff6018SBarry Smith for (i=0; i<m; i++) { 1999a0ff6018SBarry Smith row = rstart + i; 200000e6dbe6SBarry Smith nz = ii[i+1] - ii[i]; 200100e6dbe6SBarry Smith cwork = jj; jj += nz; 200200e6dbe6SBarry Smith vwork = aa; aa += nz; 20038c638d02SBarry Smith ierr = MatSetValues_MPIAIJ(M,1,&row,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 2004a0ff6018SBarry Smith } 2005a0ff6018SBarry Smith 2006a0ff6018SBarry Smith ierr = MatDestroyMatrices(1,&local);CHKERRQ(ierr); 2007a0ff6018SBarry Smith ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2008a0ff6018SBarry Smith ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2009a0ff6018SBarry Smith *newmat = M; 2010a0ff6018SBarry Smith PetscFunctionReturn(0); 2011a0ff6018SBarry Smith } 2012