1b57c8cebSBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*5cd90555SBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.234 1998/04/03 23:15:06 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; \ 57a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero into 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 \ 63a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,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; \ 128a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero into 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 \ 134a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,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) 206a8c6a408SBarry Smith if (im[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 207a8c6a408SBarry Smith if (im[i] >= aij->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,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) 219a8c6a408SBarry Smith else if (in[j] < 0) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");} 220a8c6a408SBarry Smith else if (in[j] >= aij->N) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column 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++ ) { 272a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 273a8c6a408SBarry Smith if (idxm[i] >= aij->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,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++ ) { 277a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 278a8c6a408SBarry Smith if (idxn[j] >= aij->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column 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 } 294a8c6a408SBarry Smith } else { 295a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 296b49de8d1SLois Curfman McInnes } 297b49de8d1SLois Curfman McInnes } 2983a40ed3dSBarry Smith PetscFunctionReturn(0); 299b49de8d1SLois Curfman McInnes } 300b49de8d1SLois Curfman McInnes 3015615d1e5SSatish Balay #undef __FUNC__ 3025615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ" 3038f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode) 3048a729477SBarry Smith { 30544a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 306d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 30717699dbbSLois Curfman McInnes int size = aij->size, *owners = aij->rowners; 30817699dbbSLois Curfman McInnes int rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr; 3091eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3106abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 3111eb62cbbSBarry Smith InsertMode addv; 3121eb62cbbSBarry Smith Scalar *rvalues,*svalues; 3131eb62cbbSBarry Smith 3143a40ed3dSBarry Smith PetscFunctionBegin; 3151eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 316ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 317dbb450caSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 318a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Some processors inserted others added"); 3191eb62cbbSBarry Smith } 32047794344SBarry Smith mat->insertmode = addv; /* in case this processor had no cache */ 3211eb62cbbSBarry Smith 3221eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3230452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 324cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 3250452661fSBarry Smith owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 3261eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3271eb62cbbSBarry Smith idx = aij->stash.idx[i]; 32817699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 3291eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3301eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 3318a729477SBarry Smith } 3328a729477SBarry Smith } 3338a729477SBarry Smith } 33417699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 3351eb62cbbSBarry Smith 3361eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3370452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 338ca161407SBarry Smith ierr = MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 33917699dbbSLois Curfman McInnes nreceives = work[rank]; 340ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 34117699dbbSLois Curfman McInnes nmax = work[rank]; 3420452661fSBarry Smith PetscFree(work); 3431eb62cbbSBarry Smith 3441eb62cbbSBarry Smith /* post receives: 3451eb62cbbSBarry Smith 1) each message will consist of ordered pairs 3461eb62cbbSBarry Smith (global index,value) we store the global index as a double 347d6dfbf8fSBarry Smith to simplify the message passing. 3481eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 3491eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 3501eb62cbbSBarry Smith this is a lot of wasted space. 3511eb62cbbSBarry Smith 3521eb62cbbSBarry Smith 3531eb62cbbSBarry Smith This could be done better. 3541eb62cbbSBarry Smith */ 355ca161407SBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues); 356ca161407SBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 3571eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 358ca161407SBarry Smith ierr = MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 359ca161407SBarry Smith comm,recv_waits+i);CHKERRQ(ierr); 3601eb62cbbSBarry Smith } 3611eb62cbbSBarry Smith 3621eb62cbbSBarry Smith /* do sends: 3631eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3641eb62cbbSBarry Smith the ith processor 3651eb62cbbSBarry Smith */ 3660452661fSBarry Smith svalues = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 367ca161407SBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 3680452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 3691eb62cbbSBarry Smith starts[0] = 0; 37017699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3711eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3721eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 3731eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 3741eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 3751eb62cbbSBarry Smith } 3760452661fSBarry Smith PetscFree(owner); 3771eb62cbbSBarry Smith starts[0] = 0; 37817699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3791eb62cbbSBarry Smith count = 0; 38017699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 3811eb62cbbSBarry Smith if (procs[i]) { 382ca161407SBarry Smith ierr = MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3831eb62cbbSBarry Smith } 3841eb62cbbSBarry Smith } 3850452661fSBarry Smith PetscFree(starts); PetscFree(nprocs); 3861eb62cbbSBarry Smith 3871eb62cbbSBarry Smith /* Free cache space */ 38855908cccSLois Curfman McInnes PLogInfo(mat,"MatAssemblyBegin_MPIAIJ:Number of off-processor values %d\n",aij->stash.n); 38978b31e54SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 3901eb62cbbSBarry Smith 3911eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 3921eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 3931eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 3941eb62cbbSBarry Smith aij->rmax = nmax; 3951eb62cbbSBarry Smith 3963a40ed3dSBarry Smith PetscFunctionReturn(0); 3971eb62cbbSBarry Smith } 39844a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat); 3991eb62cbbSBarry Smith 4005615d1e5SSatish Balay #undef __FUNC__ 4015615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ" 4028f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode) 4031eb62cbbSBarry Smith { 40444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 4051eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 406416022c9SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr; 407905e6a2fSBarry Smith int row,col,other_disassembled; 4081eb62cbbSBarry Smith Scalar *values,val; 40947794344SBarry Smith InsertMode addv = mat->insertmode; 4101eb62cbbSBarry Smith 4113a40ed3dSBarry Smith PetscFunctionBegin; 4121eb62cbbSBarry Smith /* wait on receives */ 4131eb62cbbSBarry Smith while (count) { 414ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4151eb62cbbSBarry Smith /* unpack receives into our local space */ 416d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 417ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_SCALAR,&n);CHKERRQ(ierr); 4181eb62cbbSBarry Smith n = n/3; 4191eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 420227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - aij->rstart; 421227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 4221eb62cbbSBarry Smith val = values[3*i+2]; 4231eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 4241eb62cbbSBarry Smith col -= aij->cstart; 4256fd7127cSSatish Balay ierr = MatSetValues(aij->A,1,&row,1,&col,&val,addv); CHKERRQ(ierr); 4263a40ed3dSBarry Smith } else { 42755a1b374SBarry Smith if (mat->was_assembled || mat->assembled) { 428905e6a2fSBarry Smith if (!aij->colmap) { 429905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat); CHKERRQ(ierr); 430905e6a2fSBarry Smith } 431905e6a2fSBarry Smith col = aij->colmap[col] - 1; 432ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 4332493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 434227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 435d6dfbf8fSBarry Smith } 4369e25ed09SBarry Smith } 4376fd7127cSSatish Balay ierr = MatSetValues(aij->B,1,&row,1,&col,&val,addv); CHKERRQ(ierr); 4381eb62cbbSBarry Smith } 4391eb62cbbSBarry Smith } 4401eb62cbbSBarry Smith count--; 4411eb62cbbSBarry Smith } 4420452661fSBarry Smith PetscFree(aij->recv_waits); PetscFree(aij->rvalues); 4431eb62cbbSBarry Smith 4441eb62cbbSBarry Smith /* wait on sends */ 4451eb62cbbSBarry Smith if (aij->nsends) { 4460a198c4cSBarry Smith send_status = (MPI_Status *) PetscMalloc(aij->nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 447ca161407SBarry Smith ierr = MPI_Waitall(aij->nsends,aij->send_waits,send_status);CHKERRQ(ierr); 4480452661fSBarry Smith PetscFree(send_status); 4491eb62cbbSBarry Smith } 4500452661fSBarry Smith PetscFree(aij->send_waits); PetscFree(aij->svalues); 4511eb62cbbSBarry Smith 45278b31e54SBarry Smith ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr); 45378b31e54SBarry Smith ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr); 4541eb62cbbSBarry Smith 4552493cbb0SBarry Smith /* determine if any processor has disassembled, if so we must 4562493cbb0SBarry Smith also disassemble ourselfs, in order that we may reassemble. */ 457ca161407SBarry Smith ierr = MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);CHKERRQ(ierr); 458227d817aSBarry Smith if (mat->was_assembled && !other_disassembled) { 4592493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 4602493cbb0SBarry Smith } 4612493cbb0SBarry Smith 4626d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 46378b31e54SBarry Smith ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr); 4645e42470aSBarry Smith } 46578b31e54SBarry Smith ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr); 46678b31e54SBarry Smith ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr); 4675e42470aSBarry Smith 4687a0afa10SBarry Smith if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;} 4693a40ed3dSBarry Smith PetscFunctionReturn(0); 4708a729477SBarry Smith } 4718a729477SBarry Smith 4725615d1e5SSatish Balay #undef __FUNC__ 4735615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ" 4748f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A) 4751eb62cbbSBarry Smith { 47644a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 477dbd7a890SLois Curfman McInnes int ierr; 4783a40ed3dSBarry Smith 4793a40ed3dSBarry Smith PetscFunctionBegin; 48078b31e54SBarry Smith ierr = MatZeroEntries(l->A); CHKERRQ(ierr); 48178b31e54SBarry Smith ierr = MatZeroEntries(l->B); CHKERRQ(ierr); 4823a40ed3dSBarry Smith PetscFunctionReturn(0); 4831eb62cbbSBarry Smith } 4841eb62cbbSBarry Smith 4855615d1e5SSatish Balay #undef __FUNC__ 4865615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ" 4878f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 4881eb62cbbSBarry Smith { 48944a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 49017699dbbSLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 4916a5c57faSSatish Balay int *procs,*nprocs,j,found,idx,nsends,*work,row; 49217699dbbSLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 4935392566eSBarry Smith int *rvalues,tag = A->tag,count,base,slen,n,*source; 4946a5c57faSSatish Balay int *lens,imdex,*lrows,*values,rstart=l->rstart; 495d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 4961eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 4971eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 4981eb62cbbSBarry Smith IS istmp; 4991eb62cbbSBarry Smith 5003a40ed3dSBarry Smith PetscFunctionBegin; 50177c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 50278b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5031eb62cbbSBarry Smith 5041eb62cbbSBarry Smith /* first count number of contributors to each processor */ 5050452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 506cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 5070452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 5081eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5091eb62cbbSBarry Smith idx = rows[i]; 5101eb62cbbSBarry Smith found = 0; 51117699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 5121eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 5131eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 5141eb62cbbSBarry Smith } 5151eb62cbbSBarry Smith } 516a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 5171eb62cbbSBarry Smith } 51817699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 5191eb62cbbSBarry Smith 5201eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 5210452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 522ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 52317699dbbSLois Curfman McInnes nrecvs = work[rank]; 524ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 52517699dbbSLois Curfman McInnes nmax = work[rank]; 5260452661fSBarry Smith PetscFree(work); 5271eb62cbbSBarry Smith 5281eb62cbbSBarry Smith /* post receives: */ 5293a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 530ca161407SBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 5311eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 532ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 5331eb62cbbSBarry Smith } 5341eb62cbbSBarry Smith 5351eb62cbbSBarry Smith /* do sends: 5361eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 5371eb62cbbSBarry Smith the ith processor 5381eb62cbbSBarry Smith */ 5390452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 5403a40ed3dSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 5410452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 5421eb62cbbSBarry Smith starts[0] = 0; 54317699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5441eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5451eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 5461eb62cbbSBarry Smith } 5471eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 5481eb62cbbSBarry Smith 5491eb62cbbSBarry Smith starts[0] = 0; 55017699dbbSLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5511eb62cbbSBarry Smith count = 0; 55217699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 5531eb62cbbSBarry Smith if (procs[i]) { 554ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 5551eb62cbbSBarry Smith } 5561eb62cbbSBarry Smith } 5570452661fSBarry Smith PetscFree(starts); 5581eb62cbbSBarry Smith 55917699dbbSLois Curfman McInnes base = owners[rank]; 5601eb62cbbSBarry Smith 5611eb62cbbSBarry Smith /* wait on receives */ 5620452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 5631eb62cbbSBarry Smith source = lens + nrecvs; 5641eb62cbbSBarry Smith count = nrecvs; slen = 0; 5651eb62cbbSBarry Smith while (count) { 566ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 5671eb62cbbSBarry Smith /* unpack receives into our local space */ 568ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 569d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 570d6dfbf8fSBarry Smith lens[imdex] = n; 5711eb62cbbSBarry Smith slen += n; 5721eb62cbbSBarry Smith count--; 5731eb62cbbSBarry Smith } 5740452661fSBarry Smith PetscFree(recv_waits); 5751eb62cbbSBarry Smith 5761eb62cbbSBarry Smith /* move the data into the send scatter */ 5770452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 5781eb62cbbSBarry Smith count = 0; 5791eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 5801eb62cbbSBarry Smith values = rvalues + i*nmax; 5811eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 5821eb62cbbSBarry Smith lrows[count++] = values[j] - base; 5831eb62cbbSBarry Smith } 5841eb62cbbSBarry Smith } 5850452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 5860452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 5871eb62cbbSBarry Smith 5881eb62cbbSBarry Smith /* actually zap the local rows */ 589029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 590464493b3SBarry Smith PLogObjectParent(A,istmp); 5916a5c57faSSatish Balay 5926a5c57faSSatish Balay ierr = MatZeroRows(l->A,istmp,0); CHKERRQ(ierr); 59378b31e54SBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr); 59478b31e54SBarry Smith ierr = ISDestroy(istmp); CHKERRQ(ierr); 5951eb62cbbSBarry Smith 5966a5c57faSSatish Balay if (diag) { 5976a5c57faSSatish Balay for ( i = 0; i < slen; i++ ) { 5986a5c57faSSatish Balay row = lrows[i] + rstart; 5996a5c57faSSatish Balay MatSetValues(A,1,&row,1,&row,diag,INSERT_VALUES); 6006a5c57faSSatish Balay } 6016a5c57faSSatish Balay MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); 6026a5c57faSSatish Balay MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); 6036a5c57faSSatish Balay } 6046a5c57faSSatish Balay PetscFree(lrows); 6051eb62cbbSBarry Smith /* wait on sends */ 6061eb62cbbSBarry Smith if (nsends) { 607ca161407SBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 608ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 6090452661fSBarry Smith PetscFree(send_status); 6101eb62cbbSBarry Smith } 6110452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 6121eb62cbbSBarry Smith 6133a40ed3dSBarry Smith PetscFunctionReturn(0); 6141eb62cbbSBarry Smith } 6151eb62cbbSBarry Smith 6165615d1e5SSatish Balay #undef __FUNC__ 6175615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ" 6188f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy) 6191eb62cbbSBarry Smith { 620416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 621fbd6ef76SBarry Smith int ierr,nt; 622416022c9SBarry Smith 6233a40ed3dSBarry Smith PetscFunctionBegin; 624a2ce50c7SBarry Smith ierr = VecGetLocalSize(xx,&nt); CHKERRQ(ierr); 625fbd6ef76SBarry Smith if (nt != a->n) { 626a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Incompatible partition of A and xx"); 627fbd6ef76SBarry Smith } 62843a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 629f830108cSBarry Smith ierr = (*a->A->ops->mult)(a->A,xx,yy); CHKERRQ(ierr); 63043a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 631f830108cSBarry Smith ierr = (*a->B->ops->multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr); 6323a40ed3dSBarry Smith PetscFunctionReturn(0); 6331eb62cbbSBarry Smith } 6341eb62cbbSBarry Smith 6355615d1e5SSatish Balay #undef __FUNC__ 6365615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ" 6378f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 638da3a660dSBarry Smith { 639416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 640da3a660dSBarry Smith int ierr; 6413a40ed3dSBarry Smith 6423a40ed3dSBarry Smith PetscFunctionBegin; 64343a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 644f830108cSBarry Smith ierr = (*a->A->ops->multadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 64543a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 646f830108cSBarry Smith ierr = (*a->B->ops->multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr); 6473a40ed3dSBarry Smith PetscFunctionReturn(0); 648da3a660dSBarry Smith } 649da3a660dSBarry Smith 6505615d1e5SSatish Balay #undef __FUNC__ 6515615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ" 6528f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy) 653da3a660dSBarry Smith { 654416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 655da3a660dSBarry Smith int ierr; 656da3a660dSBarry Smith 6573a40ed3dSBarry Smith PetscFunctionBegin; 658da3a660dSBarry Smith /* do nondiagonal part */ 659f830108cSBarry Smith ierr = (*a->B->ops->multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 660da3a660dSBarry Smith /* send it on its way */ 661537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 662da3a660dSBarry Smith /* do local part */ 663f830108cSBarry Smith ierr = (*a->A->ops->multtrans)(a->A,xx,yy); CHKERRQ(ierr); 664da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 665da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 666da3a660dSBarry Smith /* but is not perhaps always true. */ 667537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 6683a40ed3dSBarry Smith PetscFunctionReturn(0); 669da3a660dSBarry Smith } 670da3a660dSBarry Smith 6715615d1e5SSatish Balay #undef __FUNC__ 6725615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ" 6738f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 674da3a660dSBarry Smith { 675416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 676da3a660dSBarry Smith int ierr; 677da3a660dSBarry Smith 6783a40ed3dSBarry Smith PetscFunctionBegin; 679da3a660dSBarry Smith /* do nondiagonal part */ 680f830108cSBarry Smith ierr = (*a->B->ops->multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 681da3a660dSBarry Smith /* send it on its way */ 682537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 683da3a660dSBarry Smith /* do local part */ 684f830108cSBarry Smith ierr = (*a->A->ops->multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 685da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 686da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 687da3a660dSBarry Smith /* but is not perhaps always true. */ 6880a198c4cSBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 6893a40ed3dSBarry Smith PetscFunctionReturn(0); 690da3a660dSBarry Smith } 691da3a660dSBarry Smith 6921eb62cbbSBarry Smith /* 6931eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 6941eb62cbbSBarry Smith diagonal block 6951eb62cbbSBarry Smith */ 6965615d1e5SSatish Balay #undef __FUNC__ 6975615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ" 6988f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v) 6991eb62cbbSBarry Smith { 7003a40ed3dSBarry Smith int ierr; 701416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 7023a40ed3dSBarry Smith 7033a40ed3dSBarry Smith PetscFunctionBegin; 704a8c6a408SBarry Smith if (a->M != a->N) SETERRQ(PETSC_ERR_SUP,0,"Supports only square matrix where A->A is diag block"); 7055baf8537SBarry Smith if (a->rstart != a->cstart || a->rend != a->cend) { 706a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"row partition must equal col partition"); 7073a40ed3dSBarry Smith } 7083a40ed3dSBarry Smith ierr = MatGetDiagonal(a->A,v);CHKERRQ(ierr); 7093a40ed3dSBarry Smith PetscFunctionReturn(0); 7101eb62cbbSBarry Smith } 7111eb62cbbSBarry Smith 7125615d1e5SSatish Balay #undef __FUNC__ 7135615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ" 7148f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A) 715052efed2SBarry Smith { 716052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 717052efed2SBarry Smith int ierr; 7183a40ed3dSBarry Smith 7193a40ed3dSBarry Smith PetscFunctionBegin; 720052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 721052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 7223a40ed3dSBarry Smith PetscFunctionReturn(0); 723052efed2SBarry Smith } 724052efed2SBarry Smith 7255615d1e5SSatish Balay #undef __FUNC__ 726d4bb536fSBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" 727e1311b90SBarry Smith int MatDestroy_MPIAIJ(Mat mat) 7281eb62cbbSBarry Smith { 72944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 7301eb62cbbSBarry Smith int ierr; 73183e2fdc7SBarry Smith 7323a40ed3dSBarry Smith PetscFunctionBegin; 7333a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 734e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",aij->M,aij->N); 735a5a9c739SBarry Smith #endif 73683e2fdc7SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 7370452661fSBarry Smith PetscFree(aij->rowners); 73878b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 73978b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 7400452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 7410452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 7421eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 743a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 7447a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 7450452661fSBarry Smith PetscFree(aij); 746a5a9c739SBarry Smith PLogObjectDestroy(mat); 7470452661fSBarry Smith PetscHeaderDestroy(mat); 7483a40ed3dSBarry Smith PetscFunctionReturn(0); 7491eb62cbbSBarry Smith } 750ee50ffe9SBarry Smith 7515615d1e5SSatish Balay #undef __FUNC__ 752d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" 7538f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 7541eb62cbbSBarry Smith { 755416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 756416022c9SBarry Smith int ierr; 757416022c9SBarry Smith 7583a40ed3dSBarry Smith PetscFunctionBegin; 75917699dbbSLois Curfman McInnes if (aij->size == 1) { 760416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 761416022c9SBarry Smith } 762a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 7633a40ed3dSBarry Smith PetscFunctionReturn(0); 764416022c9SBarry Smith } 765416022c9SBarry Smith 7665615d1e5SSatish Balay #undef __FUNC__ 767d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab" 7688f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer) 769416022c9SBarry Smith { 77044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 771dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 772a56f8943SBarry Smith int ierr, format,shift = C->indexshift,rank; 773d636dbe3SBarry Smith FILE *fd; 77419bcc07fSBarry Smith ViewerType vtype; 775416022c9SBarry Smith 7763a40ed3dSBarry Smith PetscFunctionBegin; 77719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 77819bcc07fSBarry Smith if (vtype == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) { 77990ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 7800a198c4cSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 7814e220ebcSLois Curfman McInnes MatInfo info; 7824e220ebcSLois Curfman McInnes int flg; 783a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 78490ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 7854e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 78695e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 78777c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 78895e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 7894e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 79095e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 7914e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 7924e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->A,MAT_LOCAL,&info); 7934e220ebcSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used); 7944e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->B,MAT_LOCAL,&info); 7954e220ebcSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used); 796a56f8943SBarry Smith fflush(fd); 79777c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 798a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 7993a40ed3dSBarry Smith PetscFunctionReturn(0); 8003a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 8013a40ed3dSBarry Smith PetscFunctionReturn(0); 80208480c60SBarry Smith } 803416022c9SBarry Smith } 804416022c9SBarry Smith 80519bcc07fSBarry Smith if (vtype == DRAW_VIEWER) { 80619bcc07fSBarry Smith Draw draw; 80719bcc07fSBarry Smith PetscTruth isnull; 80819bcc07fSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 8093a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 81019bcc07fSBarry Smith } 81119bcc07fSBarry Smith 81219bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 81390ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 81477c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 815d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 81617699dbbSLois Curfman McInnes aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 8171eb62cbbSBarry Smith aij->cend); 81878b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 81978b31e54SBarry Smith ierr = MatView(aij->B,viewer); CHKERRQ(ierr); 820d13ab20cSBarry Smith fflush(fd); 82177c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 8223a40ed3dSBarry Smith } else { 823a56f8943SBarry Smith int size = aij->size; 824a56f8943SBarry Smith rank = aij->rank; 82517699dbbSLois Curfman McInnes if (size == 1) { 82678b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 8273a40ed3dSBarry Smith } else { 82895373324SBarry Smith /* assemble the entire matrix onto first processor. */ 82995373324SBarry Smith Mat A; 830ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 8312eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 83295373324SBarry Smith Scalar *a; 8332ee70a88SLois Curfman McInnes 83417699dbbSLois Curfman McInnes if (!rank) { 83555843e3eSBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr); 8363a40ed3dSBarry Smith } else { 83755843e3eSBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr); 83895373324SBarry Smith } 839464493b3SBarry Smith PLogObjectParent(mat,A); 840416022c9SBarry Smith 84195373324SBarry Smith /* copy over the A part */ 842ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 8432ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 84495373324SBarry Smith row = aij->rstart; 845dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 84695373324SBarry Smith for ( i=0; i<m; i++ ) { 847416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 84895373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 84995373324SBarry Smith } 8502ee70a88SLois Curfman McInnes aj = Aloc->j; 851dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 85295373324SBarry Smith 85395373324SBarry Smith /* copy over the B part */ 854ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 8552ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 85695373324SBarry Smith row = aij->rstart; 8570452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 858dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 85995373324SBarry Smith for ( i=0; i<m; i++ ) { 860416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 86195373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 86295373324SBarry Smith } 8630452661fSBarry Smith PetscFree(ct); 8646d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8656d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 86655843e3eSBarry Smith /* 86755843e3eSBarry Smith Everyone has to call to draw the matrix since the graphics waits are 86855843e3eSBarry Smith synchronized across all processors that share the Draw object 86955843e3eSBarry Smith */ 87055843e3eSBarry Smith if (!rank || vtype == DRAW_VIEWER) { 87178b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 87295373324SBarry Smith } 87378b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 87495373324SBarry Smith } 87595373324SBarry Smith } 8763a40ed3dSBarry Smith PetscFunctionReturn(0); 8771eb62cbbSBarry Smith } 8781eb62cbbSBarry Smith 8795615d1e5SSatish Balay #undef __FUNC__ 880d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ" 881e1311b90SBarry Smith int MatView_MPIAIJ(Mat mat,Viewer viewer) 882416022c9SBarry Smith { 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); 891*5cd90555SBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 8923a40ed3dSBarry Smith ierr = MatView_MPIAIJ_Binary(mat,viewer);CHKERRQ(ierr); 893*5cd90555SBarry Smith } else { 894*5cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 895416022c9SBarry Smith } 8963a40ed3dSBarry Smith PetscFunctionReturn(0); 897416022c9SBarry Smith } 898416022c9SBarry Smith 8991eb62cbbSBarry Smith /* 9001eb62cbbSBarry Smith This has to provide several versions. 9011eb62cbbSBarry Smith 9021eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 9031eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 904d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 9051eb62cbbSBarry Smith */ 9065615d1e5SSatish Balay #undef __FUNC__ 9075615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ" 9088f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 909dbb450caSBarry Smith double fshift,int its,Vec xx) 9108a729477SBarry Smith { 91144a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 912d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 913ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 914c16cb8f2SBarry Smith Scalar *b,*x,*xs,*ls,d,*v,sum; 9156abc6512SBarry Smith int ierr,*idx, *diag; 916416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 9178a729477SBarry Smith 9183a40ed3dSBarry Smith PetscFunctionBegin; 919d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 920dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 921dbb450caSBarry Smith ls = ls + shift; 92283e2fdc7SBarry Smith if (!A->diag) {ierr = MatMarkDiag_SeqAIJ(AA); CHKERRQ(ierr);} 923d6dfbf8fSBarry Smith diag = A->diag; 924c16cb8f2SBarry Smith if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 925da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 926f830108cSBarry Smith ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9273a40ed3dSBarry Smith PetscFunctionReturn(0); 928da3a660dSBarry Smith } 9293a40ed3dSBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 9303a40ed3dSBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 931d6dfbf8fSBarry Smith while (its--) { 932d6dfbf8fSBarry Smith /* go down through the rows */ 933d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 934d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 935dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 936dbb450caSBarry Smith v = A->a + A->i[i] + shift; 937d6dfbf8fSBarry Smith sum = b[i]; 938d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 939dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 940d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 941dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 942dbb450caSBarry Smith v = B->a + B->i[i] + shift; 943d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 94455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 945d6dfbf8fSBarry Smith } 946d6dfbf8fSBarry Smith /* come up through the rows */ 947d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 948d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 949dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 950dbb450caSBarry Smith v = A->a + A->i[i] + shift; 951d6dfbf8fSBarry Smith sum = b[i]; 952d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 953dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 954d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 955dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 956dbb450caSBarry Smith v = B->a + B->i[i] + shift; 957d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 95855a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 959d6dfbf8fSBarry Smith } 960d6dfbf8fSBarry Smith } 9613a40ed3dSBarry Smith } else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 962da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 963f830108cSBarry Smith ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9643a40ed3dSBarry Smith PetscFunctionReturn(0); 965da3a660dSBarry Smith } 9663a40ed3dSBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 9673a40ed3dSBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 968d6dfbf8fSBarry Smith while (its--) { 969d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 970d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 971dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 972dbb450caSBarry Smith v = A->a + A->i[i] + shift; 973d6dfbf8fSBarry Smith sum = b[i]; 974d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 975dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 976d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 977dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 978dbb450caSBarry Smith v = B->a + B->i[i] + shift; 979d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 98055a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 981d6dfbf8fSBarry Smith } 982d6dfbf8fSBarry Smith } 9833a40ed3dSBarry Smith } else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 984da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 985f830108cSBarry Smith ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9863a40ed3dSBarry Smith PetscFunctionReturn(0); 987da3a660dSBarry Smith } 98843a90d84SBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 98978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 99043a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 99178b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 992d6dfbf8fSBarry Smith while (its--) { 993d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 994d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 995dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 996dbb450caSBarry Smith v = A->a + A->i[i] + shift; 997d6dfbf8fSBarry Smith sum = b[i]; 998d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 999dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 1000d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 1001dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 1002dbb450caSBarry Smith v = B->a + B->i[i] + shift; 1003d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 100455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 1005d6dfbf8fSBarry Smith } 1006d6dfbf8fSBarry Smith } 10073a40ed3dSBarry Smith } else { 1008a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Parallel SOR not supported"); 1009c16cb8f2SBarry Smith } 10103a40ed3dSBarry Smith PetscFunctionReturn(0); 10118a729477SBarry Smith } 1012a66be287SLois Curfman McInnes 10135615d1e5SSatish Balay #undef __FUNC__ 1014d4bb536fSBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" 10158f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info) 1016a66be287SLois Curfman McInnes { 1017a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1018a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 10194e220ebcSLois Curfman McInnes int ierr; 10204e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 1021a66be287SLois Curfman McInnes 10223a40ed3dSBarry Smith PetscFunctionBegin; 10234e220ebcSLois Curfman McInnes info->block_size = 1.0; 10244e220ebcSLois Curfman McInnes ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr); 10254e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 10264e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 10274e220ebcSLois Curfman McInnes ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr); 10284e220ebcSLois Curfman McInnes isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded; 10294e220ebcSLois Curfman McInnes isend[3] += info->memory; isend[4] += info->mallocs; 1030a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 10314e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 10324e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 10334e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 10344e220ebcSLois Curfman McInnes info->memory = isend[3]; 10354e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 1036a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 1037ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm);CHKERRQ(ierr); 10384e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10394e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10404e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10414e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10424e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1043a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 1044ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm);CHKERRQ(ierr); 10454e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10464e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10474e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10484e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10494e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1050a66be287SLois Curfman McInnes } 10514e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 10524e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10534e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10548d700155SBarry Smith info->rows_global = (double)mat->M; 10558d700155SBarry Smith info->columns_global = (double)mat->N; 10568d700155SBarry Smith info->rows_local = (double)mat->m; 10578d700155SBarry Smith info->columns_local = (double)mat->N; 10584e220ebcSLois Curfman McInnes 10593a40ed3dSBarry Smith PetscFunctionReturn(0); 1060a66be287SLois Curfman McInnes } 1061a66be287SLois Curfman McInnes 10625615d1e5SSatish Balay #undef __FUNC__ 1063d4bb536fSBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" 10648f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op) 1065c74985f6SBarry Smith { 1066c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1067c74985f6SBarry Smith 10683a40ed3dSBarry Smith PetscFunctionBegin; 10696d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 10706d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 10716da5968aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 1072c2653b3dSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 107396854ed6SLois Curfman McInnes op == MAT_NEW_NONZERO_ALLOCATION_ERROR || 1074c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR) { 1075b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 1076b1fbbac0SLois Curfman McInnes MatSetOption(a->B,op); 1077b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 1078aeafbbfcSLois Curfman McInnes a->roworiented = 1; 1079c0bbcb79SLois Curfman McInnes MatSetOption(a->A,op); 1080c0bbcb79SLois Curfman McInnes MatSetOption(a->B,op); 1081b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 10826da5968aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 10836d4a8577SBarry Smith op == MAT_SYMMETRIC || 10846d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 10856d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 1086981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIAIJ:Option ignored\n"); 10876d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) { 10884b0e389bSBarry Smith a->roworiented = 0; 10894b0e389bSBarry Smith MatSetOption(a->A,op); 10904b0e389bSBarry Smith MatSetOption(a->B,op); 10912b362799SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 109290f02eecSBarry Smith a->donotstash = 1; 10933a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS){ 10943a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 10953a40ed3dSBarry Smith } else { 10963a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 10973a40ed3dSBarry Smith } 10983a40ed3dSBarry Smith PetscFunctionReturn(0); 1099c74985f6SBarry Smith } 1100c74985f6SBarry Smith 11015615d1e5SSatish Balay #undef __FUNC__ 1102d4bb536fSBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" 11038f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1104c74985f6SBarry Smith { 110544a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 11063a40ed3dSBarry Smith 11073a40ed3dSBarry Smith PetscFunctionBegin; 11080752156aSBarry Smith if (m) *m = mat->M; 11090752156aSBarry Smith if (n) *n = mat->N; 11103a40ed3dSBarry Smith PetscFunctionReturn(0); 1111c74985f6SBarry Smith } 1112c74985f6SBarry Smith 11135615d1e5SSatish Balay #undef __FUNC__ 1114d4bb536fSBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" 11158f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1116c74985f6SBarry Smith { 111744a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 11183a40ed3dSBarry Smith 11193a40ed3dSBarry Smith PetscFunctionBegin; 11200752156aSBarry Smith if (m) *m = mat->m; 1121f830108cSBarry Smith if (n) *n = mat->n; 11223a40ed3dSBarry Smith PetscFunctionReturn(0); 1123c74985f6SBarry Smith } 1124c74985f6SBarry Smith 11255615d1e5SSatish Balay #undef __FUNC__ 1126d4bb536fSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" 11278f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1128c74985f6SBarry Smith { 112944a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 11303a40ed3dSBarry Smith 11313a40ed3dSBarry Smith PetscFunctionBegin; 1132c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 11333a40ed3dSBarry Smith PetscFunctionReturn(0); 1134c74985f6SBarry Smith } 1135c74985f6SBarry Smith 11366d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11376d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11386d84be18SBarry Smith 11395615d1e5SSatish Balay #undef __FUNC__ 11405615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ" 11416d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 114239e00950SLois Curfman McInnes { 1143154123eaSLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 114470f0671dSBarry Smith Scalar *vworkA, *vworkB, **pvA, **pvB,*v_p; 1145154123eaSLois Curfman McInnes int i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart; 1146154123eaSLois Curfman McInnes int nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend; 114770f0671dSBarry Smith int *cmap, *idx_p; 114839e00950SLois Curfman McInnes 11493a40ed3dSBarry Smith PetscFunctionBegin; 1150a8c6a408SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Already active"); 11517a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11527a0afa10SBarry Smith 115370f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11547a0afa10SBarry Smith /* 11557a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11567a0afa10SBarry Smith */ 11577a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 1158c16cb8f2SBarry Smith int max = 1,m = mat->m,tmp; 1159c16cb8f2SBarry Smith for ( i=0; i<m; i++ ) { 11607a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11617a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11627a0afa10SBarry Smith } 11637a0afa10SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar))); 11647a0afa10SBarry Smith CHKPTRQ(mat->rowvalues); 11657a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11667a0afa10SBarry Smith } 11677a0afa10SBarry Smith 1168a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Only local rows") 1169abc0e9e4SLois Curfman McInnes lrow = row - rstart; 117039e00950SLois Curfman McInnes 1171154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1172154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1173154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 1174f830108cSBarry Smith ierr = (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 1175f830108cSBarry Smith ierr = (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1176154123eaSLois Curfman McInnes nztot = nzA + nzB; 1177154123eaSLois Curfman McInnes 117870f0671dSBarry Smith cmap = mat->garray; 1179154123eaSLois Curfman McInnes if (v || idx) { 1180154123eaSLois Curfman McInnes if (nztot) { 1181154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 118270f0671dSBarry Smith int imark = -1; 1183154123eaSLois Curfman McInnes if (v) { 118470f0671dSBarry Smith *v = v_p = mat->rowvalues; 118539e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 118670f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1187154123eaSLois Curfman McInnes else break; 1188154123eaSLois Curfman McInnes } 1189154123eaSLois Curfman McInnes imark = i; 119070f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 119170f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1192154123eaSLois Curfman McInnes } 1193154123eaSLois Curfman McInnes if (idx) { 119470f0671dSBarry Smith *idx = idx_p = mat->rowindices; 119570f0671dSBarry Smith if (imark > -1) { 119670f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 119770f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 119870f0671dSBarry Smith } 119970f0671dSBarry Smith } else { 1200154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 120170f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1202154123eaSLois Curfman McInnes else break; 1203154123eaSLois Curfman McInnes } 1204154123eaSLois Curfman McInnes imark = i; 120570f0671dSBarry Smith } 120670f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 120770f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 120839e00950SLois Curfman McInnes } 120939e00950SLois Curfman McInnes } 12101ca473b0SSatish Balay else { 12111ca473b0SSatish Balay if (idx) *idx = 0; 12121ca473b0SSatish Balay if (v) *v = 0; 12131ca473b0SSatish Balay } 1214154123eaSLois Curfman McInnes } 121539e00950SLois Curfman McInnes *nz = nztot; 1216f830108cSBarry Smith ierr = (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 1217f830108cSBarry Smith ierr = (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 12183a40ed3dSBarry Smith PetscFunctionReturn(0); 121939e00950SLois Curfman McInnes } 122039e00950SLois Curfman McInnes 12215615d1e5SSatish Balay #undef __FUNC__ 1222d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" 12236d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 122439e00950SLois Curfman McInnes { 12257a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 12263a40ed3dSBarry Smith 12273a40ed3dSBarry Smith PetscFunctionBegin; 12287a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 1229a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"MatGetRow not called"); 12307a0afa10SBarry Smith } 12317a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 12323a40ed3dSBarry Smith PetscFunctionReturn(0); 123339e00950SLois Curfman McInnes } 123439e00950SLois Curfman McInnes 12355615d1e5SSatish Balay #undef __FUNC__ 12365615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ" 12378f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1238855ac2c5SLois Curfman McInnes { 1239855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1240ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1241416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1242416022c9SBarry Smith double sum = 0.0; 124304ca555eSLois Curfman McInnes Scalar *v; 124404ca555eSLois Curfman McInnes 12453a40ed3dSBarry Smith PetscFunctionBegin; 124617699dbbSLois Curfman McInnes if (aij->size == 1) { 124714183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 124837fa93a5SLois Curfman McInnes } else { 124904ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 125004ca555eSLois Curfman McInnes v = amat->a; 125104ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 12523a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 125304ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 125404ca555eSLois Curfman McInnes #else 125504ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 125604ca555eSLois Curfman McInnes #endif 125704ca555eSLois Curfman McInnes } 125804ca555eSLois Curfman McInnes v = bmat->a; 125904ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 12603a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 126104ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 126204ca555eSLois Curfman McInnes #else 126304ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 126404ca555eSLois Curfman McInnes #endif 126504ca555eSLois Curfman McInnes } 1266ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr); 126704ca555eSLois Curfman McInnes *norm = sqrt(*norm); 12683a40ed3dSBarry Smith } else if (type == NORM_1) { /* max column norm */ 126904ca555eSLois Curfman McInnes double *tmp, *tmp2; 127004ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 1271758f045eSSatish Balay tmp = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp); 1272758f045eSSatish Balay tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2); 1273cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 127404ca555eSLois Curfman McInnes *norm = 0.0; 127504ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 127604ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1277579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 127804ca555eSLois Curfman McInnes } 127904ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 128004ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1281579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 128204ca555eSLois Curfman McInnes } 1283ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr); 128404ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 128504ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 128604ca555eSLois Curfman McInnes } 12870452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 12883a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 1289515d9167SLois Curfman McInnes double ntemp = 0.0; 129004ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1291dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 129204ca555eSLois Curfman McInnes sum = 0.0; 129304ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1294cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 129504ca555eSLois Curfman McInnes } 1296dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 129704ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1298cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 129904ca555eSLois Curfman McInnes } 1300515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 130104ca555eSLois Curfman McInnes } 1302ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);CHKERRQ(ierr); 1303ca161407SBarry Smith } else { 1304a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 130504ca555eSLois Curfman McInnes } 130637fa93a5SLois Curfman McInnes } 13073a40ed3dSBarry Smith PetscFunctionReturn(0); 1308855ac2c5SLois Curfman McInnes } 1309855ac2c5SLois Curfman McInnes 13105615d1e5SSatish Balay #undef __FUNC__ 13115615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ" 13128f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1313b7c46309SBarry Smith { 1314b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1315dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1316416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1317b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 13183a40ed3dSBarry Smith Mat B; 1319b7c46309SBarry Smith Scalar *array; 1320b7c46309SBarry Smith 13213a40ed3dSBarry Smith PetscFunctionBegin; 1322d4bb536fSBarry Smith if (matout == PETSC_NULL && M != N) { 1323a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place"); 1324d4bb536fSBarry Smith } 1325d4bb536fSBarry Smith 1326d4bb536fSBarry Smith ierr = MatCreateMPIAIJ(A->comm,a->n,a->m,N,M,0,PETSC_NULL,0,PETSC_NULL,&B);CHKERRQ(ierr); 1327b7c46309SBarry Smith 1328b7c46309SBarry Smith /* copy over the A part */ 1329ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1330b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1331b7c46309SBarry Smith row = a->rstart; 1332dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1333b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1334416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1335b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1336b7c46309SBarry Smith } 1337b7c46309SBarry Smith aj = Aloc->j; 13384af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1339b7c46309SBarry Smith 1340b7c46309SBarry Smith /* copy over the B part */ 1341ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1342b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1343b7c46309SBarry Smith row = a->rstart; 13440452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1345dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1346b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1347416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1348b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1349b7c46309SBarry Smith } 13500452661fSBarry Smith PetscFree(ct); 13516d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13526d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13533638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13540de55854SLois Curfman McInnes *matout = B; 13550de55854SLois Curfman McInnes } else { 1356f830108cSBarry Smith PetscOps *Abops; 1357f830108cSBarry Smith struct _MatOps *Aops; 1358f830108cSBarry Smith 13590de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13600452661fSBarry Smith PetscFree(a->rowners); 13610de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13620de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 13630452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 13640452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13650de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1366a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13670452661fSBarry Smith PetscFree(a); 1368f830108cSBarry Smith 1369f830108cSBarry Smith /* 1370f830108cSBarry Smith This is horrible, horrible code. We need to keep the 1371f830108cSBarry Smith A pointers for the bops and ops but copy everything 1372f830108cSBarry Smith else from C. 1373f830108cSBarry Smith */ 1374f830108cSBarry Smith Abops = A->bops; 1375f830108cSBarry Smith Aops = A->ops; 1376f09e8eb9SSatish Balay PetscMemcpy(A,B,sizeof(struct _p_Mat)); 1377f830108cSBarry Smith A->bops = Abops; 1378f830108cSBarry Smith A->ops = Aops; 13790452661fSBarry Smith PetscHeaderDestroy(B); 13800de55854SLois Curfman McInnes } 13813a40ed3dSBarry Smith PetscFunctionReturn(0); 1382b7c46309SBarry Smith } 1383b7c46309SBarry Smith 13845615d1e5SSatish Balay #undef __FUNC__ 13855615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ" 13864b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr) 1387a008b906SSatish Balay { 13884b967eb1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 13894b967eb1SSatish Balay Mat a = aij->A, b = aij->B; 1390a008b906SSatish Balay int ierr,s1,s2,s3; 1391a008b906SSatish Balay 13923a40ed3dSBarry Smith PetscFunctionBegin; 13934b967eb1SSatish Balay ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr); 13944b967eb1SSatish Balay if (rr) { 1395e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&s1);CHKERRQ(ierr); 1396a8c6a408SBarry Smith if (s1!=s3) SETERRQ(PETSC_ERR_ARG_SIZ,0,"right vector non-conforming local size"); 13974b967eb1SSatish Balay /* Overlap communication with computation. */ 139843a90d84SBarry Smith ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1399a008b906SSatish Balay } 14004b967eb1SSatish Balay if (ll) { 1401e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&s1);CHKERRQ(ierr); 1402a8c6a408SBarry Smith if (s1!=s2) SETERRQ(PETSC_ERR_ARG_SIZ,0,"left vector non-conforming local size"); 1403f830108cSBarry Smith ierr = (*b->ops->diagonalscale)(b,ll,0); CHKERRQ(ierr); 14044b967eb1SSatish Balay } 14054b967eb1SSatish Balay /* scale the diagonal block */ 1406f830108cSBarry Smith ierr = (*a->ops->diagonalscale)(a,ll,rr); CHKERRQ(ierr); 14074b967eb1SSatish Balay 14084b967eb1SSatish Balay if (rr) { 14094b967eb1SSatish Balay /* Do a scatter end and then right scale the off-diagonal block */ 141043a90d84SBarry Smith ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1411f830108cSBarry Smith ierr = (*b->ops->diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr); 14124b967eb1SSatish Balay } 14134b967eb1SSatish Balay 14143a40ed3dSBarry Smith PetscFunctionReturn(0); 1415a008b906SSatish Balay } 1416a008b906SSatish Balay 1417a008b906SSatish Balay 1418682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat); 14195615d1e5SSatish Balay #undef __FUNC__ 1420d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" 14218f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A) 1422682d7d0cSBarry Smith { 1423682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 14243a40ed3dSBarry Smith int ierr; 1425682d7d0cSBarry Smith 14263a40ed3dSBarry Smith PetscFunctionBegin; 14273a40ed3dSBarry Smith if (!a->rank) { 14283a40ed3dSBarry Smith ierr = MatPrintHelp_SeqAIJ(a->A);CHKERRQ(ierr); 14293a40ed3dSBarry Smith } 14303a40ed3dSBarry Smith PetscFunctionReturn(0); 1431682d7d0cSBarry Smith } 1432682d7d0cSBarry Smith 14335615d1e5SSatish Balay #undef __FUNC__ 1434d4bb536fSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" 14358f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs) 14365a838052SSatish Balay { 14373a40ed3dSBarry Smith PetscFunctionBegin; 14385a838052SSatish Balay *bs = 1; 14393a40ed3dSBarry Smith PetscFunctionReturn(0); 14405a838052SSatish Balay } 14415615d1e5SSatish Balay #undef __FUNC__ 1442d4bb536fSBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" 14438f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A) 1444bb5a7306SBarry Smith { 1445bb5a7306SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1446bb5a7306SBarry Smith int ierr; 14473a40ed3dSBarry Smith 14483a40ed3dSBarry Smith PetscFunctionBegin; 1449bb5a7306SBarry Smith ierr = MatSetUnfactored(a->A); CHKERRQ(ierr); 14503a40ed3dSBarry Smith PetscFunctionReturn(0); 1451bb5a7306SBarry Smith } 1452bb5a7306SBarry Smith 1453d4bb536fSBarry Smith #undef __FUNC__ 1454d4bb536fSBarry Smith #define __FUNC__ "MatEqual_MPIAIJ" 1455d4bb536fSBarry Smith int MatEqual_MPIAIJ(Mat A, Mat B, PetscTruth *flag) 1456d4bb536fSBarry Smith { 1457d4bb536fSBarry Smith Mat_MPIAIJ *matB = (Mat_MPIAIJ *) B->data,*matA = (Mat_MPIAIJ *) A->data; 1458d4bb536fSBarry Smith Mat a, b, c, d; 1459d4bb536fSBarry Smith PetscTruth flg; 1460d4bb536fSBarry Smith int ierr; 1461d4bb536fSBarry Smith 14623a40ed3dSBarry Smith PetscFunctionBegin; 1463a8c6a408SBarry Smith if (B->type != MATMPIAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 1464d4bb536fSBarry Smith a = matA->A; b = matA->B; 1465d4bb536fSBarry Smith c = matB->A; d = matB->B; 1466d4bb536fSBarry Smith 1467d4bb536fSBarry Smith ierr = MatEqual(a, c, &flg); CHKERRQ(ierr); 1468d4bb536fSBarry Smith if (flg == PETSC_TRUE) { 1469d4bb536fSBarry Smith ierr = MatEqual(b, d, &flg); CHKERRQ(ierr); 1470d4bb536fSBarry Smith } 1471ca161407SBarry Smith ierr = MPI_Allreduce(&flg, flag, 1, MPI_INT, MPI_LAND, A->comm);CHKERRQ(ierr); 14723a40ed3dSBarry Smith PetscFunctionReturn(0); 1473d4bb536fSBarry Smith } 1474d4bb536fSBarry Smith 14758f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int); 14762f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 14770a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 14780a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **); 14796a6a5d1dSBarry Smith extern int MatGetSubMatrix_MPIAIJ (Mat ,IS,IS,int,MatGetSubMatrixCall,Mat *); 148000e6dbe6SBarry Smith 14818a729477SBarry Smith /* -------------------------------------------------------------------*/ 14822ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ, 148339e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 148444a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 148544a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 148636ce4990SBarry Smith 0,0, 148736ce4990SBarry Smith 0,0, 148836ce4990SBarry Smith 0,0, 148944a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1490b7c46309SBarry Smith MatTranspose_MPIAIJ, 1491d4bb536fSBarry Smith MatGetInfo_MPIAIJ,MatEqual_MPIAIJ, 1492a008b906SSatish Balay MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ, 1493ee50ffe9SBarry Smith MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ, 14941eb62cbbSBarry Smith 0, 1495299609e3SLois Curfman McInnes MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ, 149636ce4990SBarry Smith 0,0,0,0, 1497d1710a03SLois Curfman McInnes MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ, 149836ce4990SBarry Smith 0,0, 149994a9d846SBarry Smith 0,0,MatConvertSameType_MPIAIJ,0,0, 1500b49de8d1SLois Curfman McInnes 0,0,0, 1501598137ffSSatish Balay MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0, 1502052efed2SBarry Smith MatPrintHelp_MPIAIJ, 15033b2fbd54SBarry Smith MatScale_MPIAIJ,0,0,0, 15040a198c4cSBarry Smith MatGetBlockSize_MPIAIJ,0,0,0,0, 150500e6dbe6SBarry Smith MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ, 1506ca161407SBarry Smith 0,0,MatGetSubMatrix_MPIAIJ}; 150736ce4990SBarry Smith 15088a729477SBarry Smith 15095615d1e5SSatish Balay #undef __FUNC__ 15105615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ" 15111987afe7SBarry Smith /*@C 1512ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 15133a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 15143a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 15153a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 15163a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 15178a729477SBarry Smith 15188a729477SBarry Smith Input Parameters: 15191eb62cbbSBarry Smith . comm - MPI communicator 15207d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 152192e8d321SLois Curfman McInnes This value should be the same as the local size used in creating the 152292e8d321SLois Curfman McInnes y vector for the matrix-vector product y = Ax. 15231a3896d6SBarry Smith . n - This value should be the same as the local size used in creating the 15241a3896d6SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 15251a3896d6SBarry Smith calculated if N is given) For square matrices n is almost always m. 15267d3e4905SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 152792e8d321SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1528ab693e5aSLois Curfman McInnes . d_nz - number of nonzeros per row in diagonal portion of local submatrix 1529ff756334SLois Curfman McInnes (same for all local rows) 15302bd5e0b2SLois Curfman McInnes . d_nzz - array containing the number of nonzeros in the various rows of the 153192e8d321SLois Curfman McInnes diagonal portion of the local submatrix (possibly different for each row) 15322bd5e0b2SLois Curfman McInnes or PETSC_NULL. You must leave room for the diagonal entry even if 15332bd5e0b2SLois Curfman McInnes it is zero. 15342bd5e0b2SLois Curfman McInnes . o_nz - number of nonzeros per row in the off-diagonal portion of local 1535ab693e5aSLois Curfman McInnes submatrix (same for all local rows). 15362bd5e0b2SLois Curfman McInnes . o_nzz - array containing the number of nonzeros in the various rows of the 15372bd5e0b2SLois Curfman McInnes off-diagonal portion of the local submatrix (possibly different for 15382bd5e0b2SLois Curfman McInnes each row) or PETSC_NULL. 15398a729477SBarry Smith 1540ff756334SLois Curfman McInnes Output Parameter: 154144cd7ae7SLois Curfman McInnes . A - the matrix 15428a729477SBarry Smith 1543ff756334SLois Curfman McInnes Notes: 1544ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1545ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 15460002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 15470002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1548ff756334SLois Curfman McInnes 1549ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1550ff756334SLois Curfman McInnes (possibly both). 1551ff756334SLois Curfman McInnes 15525511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 15535511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 15545511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 15555511cfe3SLois Curfman McInnes 15565511cfe3SLois Curfman McInnes Options Database Keys: 15575511cfe3SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 15586e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 15596e62573dSLois Curfman McInnes $ (max limit=5) 15606e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 15616e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 15626e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 15635511cfe3SLois Curfman McInnes 1564e0245417SLois Curfman McInnes Storage Information: 1565e0245417SLois Curfman McInnes For a square global matrix we define each processor's diagonal portion 1566e0245417SLois Curfman McInnes to be its local rows and the corresponding columns (a square submatrix); 1567e0245417SLois Curfman McInnes each processor's off-diagonal portion encompasses the remainder of the 1568e0245417SLois Curfman McInnes local matrix (a rectangular submatrix). 1569e0245417SLois Curfman McInnes 1570e0245417SLois Curfman McInnes The user can specify preallocated storage for the diagonal part of 15715ace5be8SLois Curfman McInnes the local submatrix with either d_nz or d_nnz (not both). Set 15725ace5be8SLois Curfman McInnes d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic 15735ace5be8SLois Curfman McInnes memory allocation. Likewise, specify preallocated storage for the 15745ace5be8SLois Curfman McInnes off-diagonal part of the local submatrix with o_nz or o_nnz (not both). 1575ff756334SLois Curfman McInnes 15765511cfe3SLois Curfman McInnes Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In 15775511cfe3SLois Curfman McInnes the figure below we depict these three local rows and all columns (0-11). 15782191d07cSBarry Smith 1579b810aeb4SBarry Smith $ 0 1 2 3 4 5 6 7 8 9 10 11 1580b810aeb4SBarry Smith $ ------------------- 15815511cfe3SLois Curfman McInnes $ row 3 | o o o d d d o o o o o o 15825511cfe3SLois Curfman McInnes $ row 4 | o o o d d d o o o o o o 15835511cfe3SLois Curfman McInnes $ row 5 | o o o d d d o o o o o o 15845511cfe3SLois Curfman McInnes $ ------------------- 1585b810aeb4SBarry Smith $ 1586b810aeb4SBarry Smith 15875511cfe3SLois Curfman McInnes Thus, any entries in the d locations are stored in the d (diagonal) 15885511cfe3SLois Curfman McInnes submatrix, and any entries in the o locations are stored in the 15895511cfe3SLois Curfman McInnes o (off-diagonal) submatrix. Note that the d and the o submatrices are 15905511cfe3SLois Curfman McInnes stored simply in the MATSEQAIJ format for compressed row storage. 15915511cfe3SLois Curfman McInnes 15925511cfe3SLois Curfman McInnes Now d_nz should indicate the number of nonzeros per row in the d matrix, 15935511cfe3SLois Curfman McInnes and o_nz should indicate the number of nonzeros per row in the o matrix. 15945511cfe3SLois Curfman McInnes In general, for PDE problems in which most nonzeros are near the diagonal, 15953d323bbdSBarry Smith one expects d_nz >> o_nz. For large problems you MUST preallocate memory 159692e8d321SLois Curfman McInnes or you will get TERRIBLE performance; see the users' manual chapter on 15976da5968aSLois Curfman McInnes matrices. 15983a511b96SLois Curfman McInnes 1599dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1600ff756334SLois Curfman McInnes 1601fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 16028a729477SBarry Smith @*/ 1603e1311b90SBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A) 16048a729477SBarry Smith { 160544cd7ae7SLois Curfman McInnes Mat B; 160644cd7ae7SLois Curfman McInnes Mat_MPIAIJ *b; 160736ce4990SBarry Smith int ierr, i,sum[2],work[2],size; 1608416022c9SBarry Smith 16093a40ed3dSBarry Smith PetscFunctionBegin; 16103914022bSBarry Smith MPI_Comm_size(comm,&size); 16113914022bSBarry Smith if (size == 1) { 16123914022bSBarry Smith if (M == PETSC_DECIDE) M = m; 16133914022bSBarry Smith if (N == PETSC_DECIDE) N = n; 16143914022bSBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,d_nz,d_nnz,A); CHKERRQ(ierr); 16153a40ed3dSBarry Smith PetscFunctionReturn(0); 16163914022bSBarry Smith } 16173914022bSBarry Smith 161844cd7ae7SLois Curfman McInnes *A = 0; 1619f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIAIJ,comm,MatDestroy,MatView); 162044cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 162144cd7ae7SLois Curfman McInnes B->data = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b); 162244cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_MPIAIJ)); 1623f830108cSBarry Smith PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps)); 1624e1311b90SBarry Smith B->ops->destroy = MatDestroy_MPIAIJ; 1625e1311b90SBarry Smith B->ops->view = MatView_MPIAIJ; 162644cd7ae7SLois Curfman McInnes B->factor = 0; 162744cd7ae7SLois Curfman McInnes B->assembled = PETSC_FALSE; 162890f02eecSBarry Smith B->mapping = 0; 1629d6dfbf8fSBarry Smith 163047794344SBarry Smith B->insertmode = NOT_SET_VALUES; 16319eb4d147SSatish Balay b->size = size; 163244cd7ae7SLois Curfman McInnes MPI_Comm_rank(comm,&b->rank); 16331eb62cbbSBarry Smith 16343a40ed3dSBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) { 1635a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz"); 16363a40ed3dSBarry Smith } 16371987afe7SBarry Smith 1638dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE || N == PETSC_DECIDE) { 16391eb62cbbSBarry Smith work[0] = m; work[1] = n; 1640ca161407SBarry Smith ierr = MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm );CHKERRQ(ierr); 1641dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE) M = sum[0]; 1642dbd7a890SLois Curfman McInnes if (N == PETSC_DECIDE) N = sum[1]; 16431eb62cbbSBarry Smith } 164444cd7ae7SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);} 164544cd7ae7SLois Curfman McInnes if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);} 164644cd7ae7SLois Curfman McInnes b->m = m; B->m = m; 164744cd7ae7SLois Curfman McInnes b->n = n; B->n = n; 164844cd7ae7SLois Curfman McInnes b->N = N; B->N = N; 164944cd7ae7SLois Curfman McInnes b->M = M; B->M = M; 16501eb62cbbSBarry Smith 16511eb62cbbSBarry Smith /* build local table of row and column ownerships */ 165244cd7ae7SLois Curfman McInnes b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners); 1653f09e8eb9SSatish Balay PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1654603f58a4SSatish Balay b->cowners = b->rowners + b->size + 2; 1655ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 165644cd7ae7SLois Curfman McInnes b->rowners[0] = 0; 165744cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 165844cd7ae7SLois Curfman McInnes b->rowners[i] += b->rowners[i-1]; 16598a729477SBarry Smith } 166044cd7ae7SLois Curfman McInnes b->rstart = b->rowners[b->rank]; 166144cd7ae7SLois Curfman McInnes b->rend = b->rowners[b->rank+1]; 1662ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 166344cd7ae7SLois Curfman McInnes b->cowners[0] = 0; 166444cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 166544cd7ae7SLois Curfman McInnes b->cowners[i] += b->cowners[i-1]; 16668a729477SBarry Smith } 166744cd7ae7SLois Curfman McInnes b->cstart = b->cowners[b->rank]; 166844cd7ae7SLois Curfman McInnes b->cend = b->cowners[b->rank+1]; 16698a729477SBarry Smith 16705ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 1671029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr); 167244cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->A); 16737b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 1674029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr); 167544cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->B); 16768a729477SBarry Smith 16771eb62cbbSBarry Smith /* build cache for off array entries formed */ 167844cd7ae7SLois Curfman McInnes ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr); 167990f02eecSBarry Smith b->donotstash = 0; 168044cd7ae7SLois Curfman McInnes b->colmap = 0; 168144cd7ae7SLois Curfman McInnes b->garray = 0; 168244cd7ae7SLois Curfman McInnes b->roworiented = 1; 16838a729477SBarry Smith 16841eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 168544cd7ae7SLois Curfman McInnes b->lvec = 0; 168644cd7ae7SLois Curfman McInnes b->Mvctx = 0; 16878a729477SBarry Smith 16887a0afa10SBarry Smith /* stuff for MatGetRow() */ 168944cd7ae7SLois Curfman McInnes b->rowindices = 0; 169044cd7ae7SLois Curfman McInnes b->rowvalues = 0; 169144cd7ae7SLois Curfman McInnes b->getrowactive = PETSC_FALSE; 16927a0afa10SBarry Smith 169344cd7ae7SLois Curfman McInnes *A = B; 16943a40ed3dSBarry Smith PetscFunctionReturn(0); 1695d6dfbf8fSBarry Smith } 1696c74985f6SBarry Smith 16975615d1e5SSatish Balay #undef __FUNC__ 16985615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ" 16998f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues) 1700d6dfbf8fSBarry Smith { 1701d6dfbf8fSBarry Smith Mat mat; 1702416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 1703a1b97e82SLois Curfman McInnes int ierr, len=0, flg; 1704d6dfbf8fSBarry Smith 17053a40ed3dSBarry Smith PetscFunctionBegin; 1706416022c9SBarry Smith *newmat = 0; 1707f830108cSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIAIJ,matin->comm,MatDestroy,MatView); 1708a5a9c739SBarry Smith PLogObjectCreate(mat); 17090452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1710f830108cSBarry Smith PetscMemcpy(mat->ops,&MatOps,sizeof(struct _MatOps)); 1711e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIAIJ; 1712e1311b90SBarry Smith mat->ops->view = MatView_MPIAIJ; 1713d6dfbf8fSBarry Smith mat->factor = matin->factor; 1714c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1715d6dfbf8fSBarry Smith 171644cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 171744cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 171844cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 171944cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 1720d6dfbf8fSBarry Smith 1721416022c9SBarry Smith a->rstart = oldmat->rstart; 1722416022c9SBarry Smith a->rend = oldmat->rend; 1723416022c9SBarry Smith a->cstart = oldmat->cstart; 1724416022c9SBarry Smith a->cend = oldmat->cend; 172517699dbbSLois Curfman McInnes a->size = oldmat->size; 172617699dbbSLois Curfman McInnes a->rank = oldmat->rank; 172747794344SBarry Smith mat->insertmode = NOT_SET_VALUES; 1728bcd2baecSBarry Smith a->rowvalues = 0; 1729bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1730d6dfbf8fSBarry Smith 1731603f58a4SSatish Balay a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1732f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1733603f58a4SSatish Balay a->cowners = a->rowners + a->size + 2; 1734603f58a4SSatish Balay PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int)); 1735416022c9SBarry Smith ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 17362ee70a88SLois Curfman McInnes if (oldmat->colmap) { 17370452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1738416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1739416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1740416022c9SBarry Smith } else a->colmap = 0; 17413f41c07dSBarry Smith if (oldmat->garray) { 17423f41c07dSBarry Smith len = ((Mat_SeqAIJ *) (oldmat->B->data))->n; 17433f41c07dSBarry Smith a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray); 1744464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 17453f41c07dSBarry Smith if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1746416022c9SBarry Smith } else a->garray = 0; 1747d6dfbf8fSBarry Smith 1748416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1749416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1750a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1751416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 1752416022c9SBarry Smith ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 1753416022c9SBarry Smith PLogObjectParent(mat,a->A); 1754416022c9SBarry Smith ierr = MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr); 1755416022c9SBarry Smith PLogObjectParent(mat,a->B); 17565dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 175725cdf11fSBarry Smith if (flg) { 1758682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1759682d7d0cSBarry Smith } 17608a729477SBarry Smith *newmat = mat; 17613a40ed3dSBarry Smith PetscFunctionReturn(0); 17628a729477SBarry Smith } 1763416022c9SBarry Smith 176477c4ece6SBarry Smith #include "sys.h" 1765416022c9SBarry Smith 17665615d1e5SSatish Balay #undef __FUNC__ 17675615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ" 176819bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat) 1769416022c9SBarry Smith { 1770d65a2f8fSBarry Smith Mat A; 1771d65a2f8fSBarry Smith Scalar *vals,*svals; 177219bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 1773416022c9SBarry Smith MPI_Status status; 17743a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 177517699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1776d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 177719bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 1778416022c9SBarry Smith 17793a40ed3dSBarry Smith PetscFunctionBegin; 178017699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 178117699dbbSLois Curfman McInnes if (!rank) { 178290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 17830752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr); 1784a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 1785d64ed03dSBarry Smith if (header[3] < 0) { 1786a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix in special format on disk, cannot load as MPIAIJ"); 1787d64ed03dSBarry Smith } 17886c5fab8fSBarry Smith } 17896c5fab8fSBarry Smith 1790d64ed03dSBarry Smith 1791ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 1792416022c9SBarry Smith M = header[1]; N = header[2]; 1793416022c9SBarry Smith /* determine ownership of all rows */ 179417699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 17950452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1796ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1797416022c9SBarry Smith rowners[0] = 0; 179817699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1799416022c9SBarry Smith rowners[i] += rowners[i-1]; 1800416022c9SBarry Smith } 180117699dbbSLois Curfman McInnes rstart = rowners[rank]; 180217699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1803416022c9SBarry Smith 1804416022c9SBarry Smith /* distribute row lengths to all processors */ 18050452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1806416022c9SBarry Smith offlens = ourlens + (rend-rstart); 180717699dbbSLois Curfman McInnes if (!rank) { 18080452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 18090752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 18100452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 181117699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1812ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 18130452661fSBarry Smith PetscFree(sndcounts); 18143a40ed3dSBarry Smith } else { 1815ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 1816416022c9SBarry Smith } 1817416022c9SBarry Smith 181817699dbbSLois Curfman McInnes if (!rank) { 1819416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 18200452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1821cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 182217699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1823416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1824416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1825416022c9SBarry Smith } 1826416022c9SBarry Smith } 18270452661fSBarry Smith PetscFree(rowlengths); 1828416022c9SBarry Smith 1829416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1830416022c9SBarry Smith maxnz = 0; 183117699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 18320452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1833416022c9SBarry Smith } 18340452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1835416022c9SBarry Smith 1836416022c9SBarry Smith /* read in my part of the matrix column indices */ 1837416022c9SBarry Smith nz = procsnz[0]; 18380452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 18390752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr); 1840d65a2f8fSBarry Smith 1841d65a2f8fSBarry Smith /* read in every one elses and ship off */ 184217699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1843d65a2f8fSBarry Smith nz = procsnz[i]; 18440752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 1845ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 1846d65a2f8fSBarry Smith } 18470452661fSBarry Smith PetscFree(cols); 18483a40ed3dSBarry Smith } else { 1849416022c9SBarry Smith /* determine buffer space needed for message */ 1850416022c9SBarry Smith nz = 0; 1851416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1852416022c9SBarry Smith nz += ourlens[i]; 1853416022c9SBarry Smith } 18540452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1855416022c9SBarry Smith 1856416022c9SBarry Smith /* receive message of column indices*/ 1857ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1858ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1859a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 1860416022c9SBarry Smith } 1861416022c9SBarry Smith 1862416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 1863cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 1864416022c9SBarry Smith jj = 0; 1865416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1866416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 1867d65a2f8fSBarry Smith if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1868416022c9SBarry Smith jj++; 1869416022c9SBarry Smith } 1870416022c9SBarry Smith } 1871d65a2f8fSBarry Smith 1872d65a2f8fSBarry Smith /* create our matrix */ 1873416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1874416022c9SBarry Smith ourlens[i] -= offlens[i]; 1875416022c9SBarry Smith } 1876d65a2f8fSBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 1877d65a2f8fSBarry Smith A = *newmat; 18786d4a8577SBarry Smith MatSetOption(A,MAT_COLUMNS_SORTED); 1879d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1880d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 1881d65a2f8fSBarry Smith } 1882416022c9SBarry Smith 188317699dbbSLois Curfman McInnes if (!rank) { 18840452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1885416022c9SBarry Smith 1886416022c9SBarry Smith /* read in my part of the matrix numerical values */ 1887416022c9SBarry Smith nz = procsnz[0]; 18880752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1889d65a2f8fSBarry Smith 1890d65a2f8fSBarry Smith /* insert into matrix */ 1891d65a2f8fSBarry Smith jj = rstart; 1892d65a2f8fSBarry Smith smycols = mycols; 1893d65a2f8fSBarry Smith svals = vals; 1894d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1895d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1896d65a2f8fSBarry Smith smycols += ourlens[i]; 1897d65a2f8fSBarry Smith svals += ourlens[i]; 1898d65a2f8fSBarry Smith jj++; 1899416022c9SBarry Smith } 1900416022c9SBarry Smith 1901d65a2f8fSBarry Smith /* read in other processors and ship out */ 190217699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1903416022c9SBarry Smith nz = procsnz[i]; 19040752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1905ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 1906416022c9SBarry Smith } 19070452661fSBarry Smith PetscFree(procsnz); 19083a40ed3dSBarry Smith } else { 1909d65a2f8fSBarry Smith /* receive numeric values */ 19100452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1911416022c9SBarry Smith 1912d65a2f8fSBarry Smith /* receive message of values*/ 1913ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1914ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1915a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 1916d65a2f8fSBarry Smith 1917d65a2f8fSBarry Smith /* insert into matrix */ 1918d65a2f8fSBarry Smith jj = rstart; 1919d65a2f8fSBarry Smith smycols = mycols; 1920d65a2f8fSBarry Smith svals = vals; 1921d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1922d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1923d65a2f8fSBarry Smith smycols += ourlens[i]; 1924d65a2f8fSBarry Smith svals += ourlens[i]; 1925d65a2f8fSBarry Smith jj++; 1926d65a2f8fSBarry Smith } 1927d65a2f8fSBarry Smith } 19280452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1929d65a2f8fSBarry Smith 19306d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 19316d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 19323a40ed3dSBarry Smith PetscFunctionReturn(0); 1933416022c9SBarry Smith } 1934a0ff6018SBarry Smith 193529da9460SBarry Smith #undef __FUNC__ 193629da9460SBarry Smith #define __FUNC__ "MatGetSubMatrix_MPIAIJ" 1937a0ff6018SBarry Smith /* 193829da9460SBarry Smith Not great since it makes two copies of the submatrix, first an SeqAIJ 193929da9460SBarry Smith in local and then by concatenating the local matrices the end result. 194029da9460SBarry Smith Writing it directly would be much like MatGetSubMatrices_MPIAIJ() 1941a0ff6018SBarry Smith */ 19426a6a5d1dSBarry Smith int MatGetSubMatrix_MPIAIJ(Mat mat,IS isrow,IS iscol,int csize,MatGetSubMatrixCall call,Mat *newmat) 1943a0ff6018SBarry Smith { 194400e6dbe6SBarry Smith int ierr, i, m,n,rstart,row,rend,nz,*cwork,size,rank,j; 1945a0ff6018SBarry Smith Mat *local,M; 194600e6dbe6SBarry Smith Scalar *vwork,*aa; 194700e6dbe6SBarry Smith MPI_Comm comm = mat->comm; 194800e6dbe6SBarry Smith Mat_SeqAIJ *aij; 194900e6dbe6SBarry Smith int *ii, *jj,nlocal,*dlens,*olens,dlen,olen,jend; 1950a0ff6018SBarry Smith 1951a0ff6018SBarry Smith PetscFunctionBegin; 195200e6dbe6SBarry Smith MPI_Comm_rank(comm,&rank); 195300e6dbe6SBarry Smith MPI_Comm_size(comm,&size); 195400e6dbe6SBarry Smith 1955a0ff6018SBarry Smith ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr); 1956a0ff6018SBarry Smith 1957a0ff6018SBarry Smith /* 1958a0ff6018SBarry Smith m - number of local rows 1959a0ff6018SBarry Smith n - number of columns (same on all processors) 1960a0ff6018SBarry Smith rstart - first row in new global matrix generated 1961a0ff6018SBarry Smith */ 1962a0ff6018SBarry Smith ierr = MatGetSize(*local,&m,&n);CHKERRQ(ierr); 1963a0ff6018SBarry Smith if (call == MAT_INITIAL_MATRIX) { 196400e6dbe6SBarry Smith aij = (Mat_SeqAIJ *) (*local)->data; 1965a8c6a408SBarry Smith if (aij->indexshift) SETERRQ(PETSC_ERR_SUP,1,"No support for index shifted matrix"); 196600e6dbe6SBarry Smith ii = aij->i; 196700e6dbe6SBarry Smith jj = aij->j; 196800e6dbe6SBarry Smith 1969a0ff6018SBarry Smith /* 197000e6dbe6SBarry Smith Determine the number of non-zeros in the diagonal and off-diagonal 197100e6dbe6SBarry Smith portions of the matrix in order to do correct preallocation 1972a0ff6018SBarry Smith */ 197300e6dbe6SBarry Smith 197400e6dbe6SBarry Smith /* first get start and end of "diagonal" columns */ 19756a6a5d1dSBarry Smith if (csize == PETSC_DECIDE) { 197600e6dbe6SBarry Smith nlocal = n/size + ((n % size) > rank); 19776a6a5d1dSBarry Smith } else { 19786a6a5d1dSBarry Smith nlocal = csize; 19796a6a5d1dSBarry Smith } 1980ca161407SBarry Smith ierr = MPI_Scan(&nlocal,&rend,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 198100e6dbe6SBarry Smith rstart = rend - nlocal; 19826a6a5d1dSBarry Smith if (rank == size - 1 && rend != n) { 19836a6a5d1dSBarry Smith SETERRQ(1,1,"Local column sizes do not add up to total number of columns"); 19846a6a5d1dSBarry Smith } 198500e6dbe6SBarry Smith 198600e6dbe6SBarry Smith /* next, compute all the lengths */ 198700e6dbe6SBarry Smith dlens = (int *) PetscMalloc( (2*m+1)*sizeof(int) );CHKPTRQ(dlens); 198800e6dbe6SBarry Smith olens = dlens + m; 198900e6dbe6SBarry Smith for ( i=0; i<m; i++ ) { 199000e6dbe6SBarry Smith jend = ii[i+1] - ii[i]; 199100e6dbe6SBarry Smith olen = 0; 199200e6dbe6SBarry Smith dlen = 0; 199300e6dbe6SBarry Smith for ( j=0; j<jend; j++ ) { 199400e6dbe6SBarry Smith if ( *jj < rstart || *jj >= rend) olen++; 199500e6dbe6SBarry Smith else dlen++; 199600e6dbe6SBarry Smith jj++; 199700e6dbe6SBarry Smith } 199800e6dbe6SBarry Smith olens[i] = olen; 199900e6dbe6SBarry Smith dlens[i] = dlen; 200000e6dbe6SBarry Smith } 200100e6dbe6SBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,PETSC_DECIDE,n,0,dlens,0,olens,&M);CHKERRQ(ierr); 200200e6dbe6SBarry Smith PetscFree(dlens); 2003a0ff6018SBarry Smith } else { 2004a0ff6018SBarry Smith int ml,nl; 2005a0ff6018SBarry Smith 2006a0ff6018SBarry Smith M = *newmat; 2007a0ff6018SBarry Smith ierr = MatGetLocalSize(M,&ml,&nl);CHKERRQ(ierr); 2008a8c6a408SBarry Smith if (ml != m) SETERRQ(PETSC_ERR_ARG_SIZ,1,"Previous matrix must be same size/layout as request"); 2009a0ff6018SBarry Smith ierr = MatZeroEntries(M);CHKERRQ(ierr); 2010a0ff6018SBarry Smith } 2011a0ff6018SBarry Smith ierr = MatGetOwnershipRange(M,&rstart,&rend); CHKERRQ(ierr); 201200e6dbe6SBarry Smith aij = (Mat_SeqAIJ *) (*local)->data; 2013a8c6a408SBarry Smith if (aij->indexshift) SETERRQ(PETSC_ERR_SUP,1,"No support for index shifted matrix"); 201400e6dbe6SBarry Smith ii = aij->i; 201500e6dbe6SBarry Smith jj = aij->j; 201600e6dbe6SBarry Smith aa = aij->a; 2017a0ff6018SBarry Smith for (i=0; i<m; i++) { 2018a0ff6018SBarry Smith row = rstart + i; 201900e6dbe6SBarry Smith nz = ii[i+1] - ii[i]; 202000e6dbe6SBarry Smith cwork = jj; jj += nz; 202100e6dbe6SBarry Smith vwork = aa; aa += nz; 20228c638d02SBarry Smith ierr = MatSetValues_MPIAIJ(M,1,&row,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 2023a0ff6018SBarry Smith } 2024a0ff6018SBarry Smith 2025a0ff6018SBarry Smith ierr = MatDestroyMatrices(1,&local);CHKERRQ(ierr); 2026a0ff6018SBarry Smith ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2027a0ff6018SBarry Smith ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2028a0ff6018SBarry Smith *newmat = M; 2029a0ff6018SBarry Smith PetscFunctionReturn(0); 2030a0ff6018SBarry Smith } 2031