1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*9eb4d147SSatish Balay static char vcid[] = "$Id: mpiaij.c,v 1.212 1997/08/07 14:39:24 bsmith Exp balay $"; 3cb512458SBarry Smith #endif 48a729477SBarry Smith 53369ce9aSBarry Smith #include "pinclude/pviewer.h" 670f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 7f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 8d9942c19SSatish Balay #include "src/inline/spops.h" 98a729477SBarry Smith 109e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column 119e25ed09SBarry Smith number to the local number in the off-diagonal part of the local 129e25ed09SBarry Smith storage of the matrix. This is done in a non scable way since the 139e25ed09SBarry Smith length of colmap equals the global matrix length. 149e25ed09SBarry Smith */ 155615d1e5SSatish Balay #undef __FUNC__ 165eea60f9SBarry Smith #define __FUNC__ "CreateColmap_MPIAIJ_Private" /* ADIC Ignore */ 170a198c4cSBarry Smith int CreateColmap_MPIAIJ_Private(Mat mat) 189e25ed09SBarry Smith { 1944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 20ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data; 21905e6a2fSBarry Smith int n = B->n,i; 22dbb450caSBarry Smith 23758f045eSSatish Balay aij->colmap = (int *) PetscMalloc((aij->N+1)*sizeof(int));CHKPTRQ(aij->colmap); 24464493b3SBarry Smith PLogObjectMemory(mat,aij->N*sizeof(int)); 25cddf8d76SBarry Smith PetscMemzero(aij->colmap,aij->N*sizeof(int)); 26905e6a2fSBarry Smith for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1; 279e25ed09SBarry Smith return 0; 289e25ed09SBarry Smith } 299e25ed09SBarry Smith 302493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat); 312493cbb0SBarry Smith 320520107fSSatish Balay #define CHUNKSIZE 15 3330770e4dSSatish Balay #define MatSetValues_SeqAIJ_A_Private(row,col,value,addv) \ 340520107fSSatish Balay { \ 350520107fSSatish Balay \ 360520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 3730770e4dSSatish Balay rmax = aimax[row]; nrow = ailen[row]; \ 38f5e9677aSSatish Balay col1 = col - shift; \ 39f5e9677aSSatish Balay \ 40ba4e3ef2SSatish Balay low = 0; high = nrow; \ 41ba4e3ef2SSatish Balay while (high-low > 5) { \ 42ba4e3ef2SSatish Balay t = (low+high)/2; \ 43ba4e3ef2SSatish Balay if (rp[t] > col) high = t; \ 44ba4e3ef2SSatish Balay else low = t; \ 45ba4e3ef2SSatish Balay } \ 460520107fSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 47f5e9677aSSatish Balay if (rp[_i] > col1) break; \ 48f5e9677aSSatish Balay if (rp[_i] == col1) { \ 490520107fSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 500520107fSSatish Balay else ap[_i] = value; \ 5130770e4dSSatish Balay goto a_noinsert; \ 520520107fSSatish Balay } \ 530520107fSSatish Balay } \ 5489280ab3SLois Curfman McInnes if (nonew == 1) goto a_noinsert; \ 5511ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 560520107fSSatish Balay if (nrow >= rmax) { \ 570520107fSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 580520107fSSatish Balay int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \ 590520107fSSatish Balay Scalar *new_a; \ 600520107fSSatish Balay \ 6111ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 6289280ab3SLois Curfman McInnes \ 630520107fSSatish Balay /* malloc new storage space */ \ 640520107fSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \ 650520107fSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 660520107fSSatish Balay new_j = (int *) (new_a + new_nz); \ 670520107fSSatish Balay new_i = new_j + new_nz; \ 680520107fSSatish Balay \ 690520107fSSatish Balay /* copy over old data into new slots */ \ 700520107fSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \ 710520107fSSatish Balay for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \ 720520107fSSatish Balay PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \ 730520107fSSatish Balay len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \ 740520107fSSatish Balay PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \ 750520107fSSatish Balay len*sizeof(int)); \ 760520107fSSatish Balay PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \ 770520107fSSatish Balay PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \ 780520107fSSatish Balay len*sizeof(Scalar)); \ 790520107fSSatish Balay /* free up old matrix storage */ \ 80f5e9677aSSatish Balay \ 810520107fSSatish Balay PetscFree(a->a); \ 820520107fSSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \ 830520107fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; \ 840520107fSSatish Balay a->singlemalloc = 1; \ 850520107fSSatish Balay \ 860520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 8730770e4dSSatish Balay rmax = aimax[row] = aimax[row] + CHUNKSIZE; \ 880520107fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 890520107fSSatish Balay a->maxnz += CHUNKSIZE; \ 900520107fSSatish Balay a->reallocs++; \ 910520107fSSatish Balay } \ 920520107fSSatish Balay N = nrow++ - 1; a->nz++; \ 930520107fSSatish Balay /* shift up all the later entries in this row */ \ 940520107fSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 950520107fSSatish Balay rp[ii+1] = rp[ii]; \ 960520107fSSatish Balay ap[ii+1] = ap[ii]; \ 970520107fSSatish Balay } \ 98f5e9677aSSatish Balay rp[_i] = col1; \ 990520107fSSatish Balay ap[_i] = value; \ 10030770e4dSSatish Balay a_noinsert: ; \ 1010520107fSSatish Balay ailen[row] = nrow; \ 1020520107fSSatish Balay } 1030a198c4cSBarry Smith 10430770e4dSSatish Balay #define MatSetValues_SeqAIJ_B_Private(row,col,value,addv) \ 10530770e4dSSatish Balay { \ 10630770e4dSSatish Balay \ 10730770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 10830770e4dSSatish Balay rmax = bimax[row]; nrow = bilen[row]; \ 10930770e4dSSatish Balay col1 = col - shift; \ 11030770e4dSSatish Balay \ 111ba4e3ef2SSatish Balay low = 0; high = nrow; \ 112ba4e3ef2SSatish Balay while (high-low > 5) { \ 113ba4e3ef2SSatish Balay t = (low+high)/2; \ 114ba4e3ef2SSatish Balay if (rp[t] > col) high = t; \ 115ba4e3ef2SSatish Balay else low = t; \ 116ba4e3ef2SSatish Balay } \ 11730770e4dSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 11830770e4dSSatish Balay if (rp[_i] > col1) break; \ 11930770e4dSSatish Balay if (rp[_i] == col1) { \ 12030770e4dSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 12130770e4dSSatish Balay else ap[_i] = value; \ 12230770e4dSSatish Balay goto b_noinsert; \ 12330770e4dSSatish Balay } \ 12430770e4dSSatish Balay } \ 12589280ab3SLois Curfman McInnes if (nonew == 1) goto b_noinsert; \ 12611ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 12730770e4dSSatish Balay if (nrow >= rmax) { \ 12830770e4dSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 12974c639caSSatish Balay int new_nz = bi[b->m] + CHUNKSIZE,len,*new_i,*new_j; \ 13030770e4dSSatish Balay Scalar *new_a; \ 13130770e4dSSatish Balay \ 13211ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 13389280ab3SLois Curfman McInnes \ 13430770e4dSSatish Balay /* malloc new storage space */ \ 13574c639caSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(b->m+1)*sizeof(int); \ 13630770e4dSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 13730770e4dSSatish Balay new_j = (int *) (new_a + new_nz); \ 13830770e4dSSatish Balay new_i = new_j + new_nz; \ 13930770e4dSSatish Balay \ 14030770e4dSSatish Balay /* copy over old data into new slots */ \ 14130770e4dSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = bi[ii];} \ 14274c639caSSatish Balay for ( ii=row+1; ii<b->m+1; ii++ ) {new_i[ii] = bi[ii]+CHUNKSIZE;} \ 14330770e4dSSatish Balay PetscMemcpy(new_j,bj,(bi[row]+nrow+shift)*sizeof(int)); \ 14430770e4dSSatish Balay len = (new_nz - CHUNKSIZE - bi[row] - nrow - shift); \ 14530770e4dSSatish Balay PetscMemcpy(new_j+bi[row]+shift+nrow+CHUNKSIZE,bj+bi[row]+shift+nrow, \ 14630770e4dSSatish Balay len*sizeof(int)); \ 14730770e4dSSatish Balay PetscMemcpy(new_a,ba,(bi[row]+nrow+shift)*sizeof(Scalar)); \ 14830770e4dSSatish Balay PetscMemcpy(new_a+bi[row]+shift+nrow+CHUNKSIZE,ba+bi[row]+shift+nrow, \ 14930770e4dSSatish Balay len*sizeof(Scalar)); \ 15030770e4dSSatish Balay /* free up old matrix storage */ \ 15130770e4dSSatish Balay \ 15274c639caSSatish Balay PetscFree(b->a); \ 15374c639caSSatish Balay if (!b->singlemalloc) {PetscFree(b->i);PetscFree(b->j);} \ 15474c639caSSatish Balay ba = b->a = new_a; bi = b->i = new_i; bj = b->j = new_j; \ 15574c639caSSatish Balay b->singlemalloc = 1; \ 15630770e4dSSatish Balay \ 15730770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 15830770e4dSSatish Balay rmax = bimax[row] = bimax[row] + CHUNKSIZE; \ 15974c639caSSatish Balay PLogObjectMemory(B,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 16074c639caSSatish Balay b->maxnz += CHUNKSIZE; \ 16174c639caSSatish Balay b->reallocs++; \ 16230770e4dSSatish Balay } \ 16374c639caSSatish Balay N = nrow++ - 1; b->nz++; \ 16430770e4dSSatish Balay /* shift up all the later entries in this row */ \ 16530770e4dSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 16630770e4dSSatish Balay rp[ii+1] = rp[ii]; \ 16730770e4dSSatish Balay ap[ii+1] = ap[ii]; \ 16830770e4dSSatish Balay } \ 16930770e4dSSatish Balay rp[_i] = col1; \ 17030770e4dSSatish Balay ap[_i] = value; \ 17130770e4dSSatish Balay b_noinsert: ; \ 17230770e4dSSatish Balay bilen[row] = nrow; \ 17330770e4dSSatish Balay } 17430770e4dSSatish Balay 1750520107fSSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode); 1765615d1e5SSatish Balay #undef __FUNC__ 1775615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIAIJ" 1788f6be9afSLois Curfman McInnes int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 1798a729477SBarry Smith { 18044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1814b0e389bSBarry Smith Scalar value; 1821eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 1831eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 184905e6a2fSBarry Smith int roworiented = aij->roworiented; 1858a729477SBarry Smith 1860520107fSSatish Balay /* Some Variables required in the macro */ 1874ee7247eSSatish Balay Mat A = aij->A; 1884ee7247eSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 18930770e4dSSatish Balay int *aimax = a->imax, *ai = a->i, *ailen = a->ilen,*aj = a->j; 19030770e4dSSatish Balay Scalar *aa = a->a; 19130770e4dSSatish Balay 19230770e4dSSatish Balay Mat B = aij->B; 19330770e4dSSatish Balay Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data; 19430770e4dSSatish Balay int *bimax = b->imax, *bi = b->i, *bilen = b->ilen,*bj = b->j; 19530770e4dSSatish Balay Scalar *ba = b->a; 19630770e4dSSatish Balay 197ba4e3ef2SSatish Balay int *rp,ii,nrow,_i,rmax, N, col1,low,high,t; 19830770e4dSSatish Balay int nonew = a->nonew,shift = a->indexshift; 19930770e4dSSatish Balay Scalar *ap; 2004ee7247eSSatish Balay 2018a729477SBarry Smith for ( i=0; i<m; i++ ) { 2020a198c4cSBarry Smith #if defined(PETSC_BOPT_g) 203e3372554SBarry Smith if (im[i] < 0) SETERRQ(1,0,"Negative row"); 204e3372554SBarry Smith if (im[i] >= aij->M) SETERRQ(1,0,"Row too large"); 2050a198c4cSBarry Smith #endif 2064b0e389bSBarry Smith if (im[i] >= rstart && im[i] < rend) { 2074b0e389bSBarry Smith row = im[i] - rstart; 2081eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 2094b0e389bSBarry Smith if (in[j] >= cstart && in[j] < cend){ 2104b0e389bSBarry Smith col = in[j] - cstart; 2114b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 21230770e4dSSatish Balay MatSetValues_SeqAIJ_A_Private(row,col,value,addv); 2130520107fSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2141eb62cbbSBarry Smith } 2150a198c4cSBarry Smith #if defined(PETSC_BOPT_g) 216e3372554SBarry Smith else if (in[j] < 0) {SETERRQ(1,0,"Negative column");} 217e3372554SBarry Smith else if (in[j] >= aij->N) {SETERRQ(1,0,"Col too large");} 2180a198c4cSBarry Smith #endif 2191eb62cbbSBarry Smith else { 220227d817aSBarry Smith if (mat->was_assembled) { 221905e6a2fSBarry Smith if (!aij->colmap) { 222905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 223905e6a2fSBarry Smith } 224905e6a2fSBarry Smith col = aij->colmap[in[j]] - 1; 225ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 2262493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 2274b0e389bSBarry Smith col = in[j]; 2289bf004c3SSatish Balay /* Reinitialize the variables required by MatSetValues_SeqAIJ_B_Private() */ 229f9508a3cSSatish Balay B = aij->B; 230f9508a3cSSatish Balay b = (Mat_SeqAIJ *) B->data; 231f9508a3cSSatish Balay bimax = b->imax; bi = b->i; bilen = b->ilen; bj = b->j; 232f9508a3cSSatish Balay ba = b->a; 233d6dfbf8fSBarry Smith } 2349e25ed09SBarry Smith } 2354b0e389bSBarry Smith else col = in[j]; 2364b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 23730770e4dSSatish Balay MatSetValues_SeqAIJ_B_Private(row,col,value,addv); 23830770e4dSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2391eb62cbbSBarry Smith } 2401eb62cbbSBarry Smith } 2411eb62cbbSBarry Smith } 2421eb62cbbSBarry Smith else { 24390f02eecSBarry Smith if (roworiented && !aij->donotstash) { 2444b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr); 2454b0e389bSBarry Smith } 2464b0e389bSBarry Smith else { 24790f02eecSBarry Smith if (!aij->donotstash) { 2484b0e389bSBarry Smith row = im[i]; 2494b0e389bSBarry Smith for ( j=0; j<n; j++ ) { 2504b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr); 2514b0e389bSBarry Smith } 2524b0e389bSBarry Smith } 2531eb62cbbSBarry Smith } 2548a729477SBarry Smith } 25590f02eecSBarry Smith } 2568a729477SBarry Smith return 0; 2578a729477SBarry Smith } 2588a729477SBarry Smith 2595615d1e5SSatish Balay #undef __FUNC__ 2605615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIAIJ" 2618f6be9afSLois Curfman McInnes int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 262b49de8d1SLois Curfman McInnes { 263b49de8d1SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 264b49de8d1SLois Curfman McInnes int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 265b49de8d1SLois Curfman McInnes int cstart = aij->cstart, cend = aij->cend,row,col; 266b49de8d1SLois Curfman McInnes 267b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 268e3372554SBarry Smith if (idxm[i] < 0) SETERRQ(1,0,"Negative row"); 269e3372554SBarry Smith if (idxm[i] >= aij->M) SETERRQ(1,0,"Row too large"); 270b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 271b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 272b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 273e3372554SBarry Smith if (idxn[j] < 0) SETERRQ(1,0,"Negative column"); 274e3372554SBarry Smith if (idxn[j] >= aij->N) SETERRQ(1,0,"Col too large"); 275b49de8d1SLois Curfman McInnes if (idxn[j] >= cstart && idxn[j] < cend){ 276b49de8d1SLois Curfman McInnes col = idxn[j] - cstart; 277b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 278b49de8d1SLois Curfman McInnes } 279b49de8d1SLois Curfman McInnes else { 280905e6a2fSBarry Smith if (!aij->colmap) { 281905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 282905e6a2fSBarry Smith } 283905e6a2fSBarry Smith col = aij->colmap[idxn[j]] - 1; 284e60e1c95SSatish Balay if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0; 285d9d09a02SSatish Balay else { 286b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 287b49de8d1SLois Curfman McInnes } 288b49de8d1SLois Curfman McInnes } 289b49de8d1SLois Curfman McInnes } 290d9d09a02SSatish Balay } 291b49de8d1SLois Curfman McInnes else { 292e3372554SBarry Smith SETERRQ(1,0,"Only local values currently supported"); 293b49de8d1SLois Curfman McInnes } 294b49de8d1SLois Curfman McInnes } 295b49de8d1SLois Curfman McInnes return 0; 296b49de8d1SLois Curfman McInnes } 297b49de8d1SLois Curfman McInnes 2985615d1e5SSatish Balay #undef __FUNC__ 2995615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ" 3008f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode) 3018a729477SBarry Smith { 30244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 303d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 30417699dbbSLois Curfman McInnes int size = aij->size, *owners = aij->rowners; 30517699dbbSLois Curfman McInnes int rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr; 3061eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3076abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 3081eb62cbbSBarry Smith InsertMode addv; 3091eb62cbbSBarry Smith Scalar *rvalues,*svalues; 3101eb62cbbSBarry Smith 3111eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 31247794344SBarry Smith MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm); 313dbb450caSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 314e3372554SBarry Smith SETERRQ(1,0,"Some processors inserted others added"); 3151eb62cbbSBarry Smith } 31647794344SBarry Smith mat->insertmode = addv; /* in case this processor had no cache */ 3171eb62cbbSBarry Smith 3181eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3190452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 320cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 3210452661fSBarry Smith owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 3221eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3231eb62cbbSBarry Smith idx = aij->stash.idx[i]; 32417699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 3251eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3261eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 3278a729477SBarry Smith } 3288a729477SBarry Smith } 3298a729477SBarry Smith } 33017699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 3311eb62cbbSBarry Smith 3321eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3330452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 33417699dbbSLois Curfman McInnes MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm); 33517699dbbSLois Curfman McInnes nreceives = work[rank]; 33617699dbbSLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 33717699dbbSLois Curfman McInnes nmax = work[rank]; 3380452661fSBarry Smith PetscFree(work); 3391eb62cbbSBarry Smith 3401eb62cbbSBarry Smith /* post receives: 3411eb62cbbSBarry Smith 1) each message will consist of ordered pairs 3421eb62cbbSBarry Smith (global index,value) we store the global index as a double 343d6dfbf8fSBarry Smith to simplify the message passing. 3441eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 3451eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 3461eb62cbbSBarry Smith this is a lot of wasted space. 3471eb62cbbSBarry Smith 3481eb62cbbSBarry Smith 3491eb62cbbSBarry Smith This could be done better. 3501eb62cbbSBarry Smith */ 3510452661fSBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 35278b31e54SBarry Smith CHKPTRQ(rvalues); 3530452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request)); 35478b31e54SBarry Smith CHKPTRQ(recv_waits); 3551eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 356416022c9SBarry Smith MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 3571eb62cbbSBarry Smith comm,recv_waits+i); 3581eb62cbbSBarry Smith } 3591eb62cbbSBarry Smith 3601eb62cbbSBarry Smith /* do sends: 3611eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3621eb62cbbSBarry Smith the ith processor 3631eb62cbbSBarry Smith */ 3640452661fSBarry Smith svalues = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 3650452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 36678b31e54SBarry Smith CHKPTRQ(send_waits); 3670452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 3681eb62cbbSBarry Smith starts[0] = 0; 36917699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3701eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3711eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 3721eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 3731eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 3741eb62cbbSBarry Smith } 3750452661fSBarry Smith PetscFree(owner); 3761eb62cbbSBarry Smith starts[0] = 0; 37717699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3781eb62cbbSBarry Smith count = 0; 37917699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 3801eb62cbbSBarry Smith if (procs[i]) { 381416022c9SBarry Smith MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag, 3821eb62cbbSBarry Smith comm,send_waits+count++); 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 3961eb62cbbSBarry Smith return 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 4111eb62cbbSBarry Smith /* wait on receives */ 4121eb62cbbSBarry Smith while (count) { 413d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status); 4141eb62cbbSBarry Smith /* unpack receives into our local space */ 415d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 416c60448a5SBarry Smith MPI_Get_count(&recv_status,MPIU_SCALAR,&n); 4171eb62cbbSBarry Smith n = n/3; 4181eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 419227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - aij->rstart; 420227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 4211eb62cbbSBarry Smith val = values[3*i+2]; 4221eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 4231eb62cbbSBarry Smith col -= aij->cstart; 4246fd7127cSSatish Balay ierr = MatSetValues(aij->A,1,&row,1,&col,&val,addv); CHKERRQ(ierr); 4251eb62cbbSBarry Smith } 4261eb62cbbSBarry 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); 4471eb62cbbSBarry Smith MPI_Waitall(aij->nsends,aij->send_waits,send_status); 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. */ 457227d817aSBarry Smith MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm); 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;} 4698a729477SBarry Smith return 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; 47878b31e54SBarry Smith ierr = MatZeroEntries(l->A); CHKERRQ(ierr); 47978b31e54SBarry Smith ierr = MatZeroEntries(l->B); CHKERRQ(ierr); 4801eb62cbbSBarry Smith return 0; 4811eb62cbbSBarry Smith } 4821eb62cbbSBarry Smith 4831eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 4841eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 4851eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 4861eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 4871eb62cbbSBarry Smith routine. 4881eb62cbbSBarry Smith */ 4895615d1e5SSatish Balay #undef __FUNC__ 4905615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ" 4918f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 4921eb62cbbSBarry Smith { 49344a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 49417699dbbSLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 4956abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 49617699dbbSLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 4975392566eSBarry Smith int *rvalues,tag = A->tag,count,base,slen,n,*source; 498d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 499d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 5001eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 5011eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 5021eb62cbbSBarry Smith IS istmp; 5031eb62cbbSBarry Smith 50477c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 50578b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5061eb62cbbSBarry Smith 5071eb62cbbSBarry Smith /* first count number of contributors to each processor */ 5080452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 509cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 5100452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 5111eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5121eb62cbbSBarry Smith idx = rows[i]; 5131eb62cbbSBarry Smith found = 0; 51417699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 5151eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 5161eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 5171eb62cbbSBarry Smith } 5181eb62cbbSBarry Smith } 519e3372554SBarry Smith if (!found) SETERRQ(1,0,"Index out of range"); 5201eb62cbbSBarry Smith } 52117699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 5221eb62cbbSBarry Smith 5231eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 5240452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 52517699dbbSLois Curfman McInnes MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm); 52617699dbbSLois Curfman McInnes nrecvs = work[rank]; 52717699dbbSLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 52817699dbbSLois Curfman McInnes nmax = work[rank]; 5290452661fSBarry Smith PetscFree(work); 5301eb62cbbSBarry Smith 5311eb62cbbSBarry Smith /* post receives: */ 5320452661fSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 53378b31e54SBarry Smith CHKPTRQ(rvalues); 5340452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request)); 53578b31e54SBarry Smith CHKPTRQ(recv_waits); 5361eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 537416022c9SBarry Smith MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i); 5381eb62cbbSBarry Smith } 5391eb62cbbSBarry Smith 5401eb62cbbSBarry Smith /* do sends: 5411eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 5421eb62cbbSBarry Smith the ith processor 5431eb62cbbSBarry Smith */ 5440452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 5450452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 54678b31e54SBarry Smith CHKPTRQ(send_waits); 5470452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 5481eb62cbbSBarry Smith starts[0] = 0; 54917699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5501eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5511eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 5521eb62cbbSBarry Smith } 5531eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 5541eb62cbbSBarry Smith 5551eb62cbbSBarry Smith starts[0] = 0; 55617699dbbSLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5571eb62cbbSBarry Smith count = 0; 55817699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 5591eb62cbbSBarry Smith if (procs[i]) { 560416022c9SBarry Smith MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++); 5611eb62cbbSBarry Smith } 5621eb62cbbSBarry Smith } 5630452661fSBarry Smith PetscFree(starts); 5641eb62cbbSBarry Smith 56517699dbbSLois Curfman McInnes base = owners[rank]; 5661eb62cbbSBarry Smith 5671eb62cbbSBarry Smith /* wait on receives */ 5680452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 5691eb62cbbSBarry Smith source = lens + nrecvs; 5701eb62cbbSBarry Smith count = nrecvs; slen = 0; 5711eb62cbbSBarry Smith while (count) { 572d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 5731eb62cbbSBarry Smith /* unpack receives into our local space */ 5741eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_INT,&n); 575d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 576d6dfbf8fSBarry Smith lens[imdex] = n; 5771eb62cbbSBarry Smith slen += n; 5781eb62cbbSBarry Smith count--; 5791eb62cbbSBarry Smith } 5800452661fSBarry Smith PetscFree(recv_waits); 5811eb62cbbSBarry Smith 5821eb62cbbSBarry Smith /* move the data into the send scatter */ 5830452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 5841eb62cbbSBarry Smith count = 0; 5851eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 5861eb62cbbSBarry Smith values = rvalues + i*nmax; 5871eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 5881eb62cbbSBarry Smith lrows[count++] = values[j] - base; 5891eb62cbbSBarry Smith } 5901eb62cbbSBarry Smith } 5910452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 5920452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 5931eb62cbbSBarry Smith 5941eb62cbbSBarry Smith /* actually zap the local rows */ 595029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 596464493b3SBarry Smith PLogObjectParent(A,istmp); 5970452661fSBarry Smith PetscFree(lrows); 59878b31e54SBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 59978b31e54SBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr); 60078b31e54SBarry Smith ierr = ISDestroy(istmp); CHKERRQ(ierr); 6011eb62cbbSBarry Smith 6021eb62cbbSBarry Smith /* wait on sends */ 6031eb62cbbSBarry Smith if (nsends) { 6040452661fSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status)); 60578b31e54SBarry Smith CHKPTRQ(send_status); 6061eb62cbbSBarry Smith MPI_Waitall(nsends,send_waits,send_status); 6070452661fSBarry Smith PetscFree(send_status); 6081eb62cbbSBarry Smith } 6090452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 6101eb62cbbSBarry Smith 6111eb62cbbSBarry Smith return 0; 6121eb62cbbSBarry Smith } 6131eb62cbbSBarry Smith 6145615d1e5SSatish Balay #undef __FUNC__ 6155615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ" 6168f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy) 6171eb62cbbSBarry Smith { 618416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 619fbd6ef76SBarry Smith int ierr,nt; 620416022c9SBarry Smith 621a2ce50c7SBarry Smith ierr = VecGetLocalSize(xx,&nt); CHKERRQ(ierr); 622fbd6ef76SBarry Smith if (nt != a->n) { 623f40265daSLois Curfman McInnes SETERRQ(1,0,"Incompatible partition of A and xx"); 624fbd6ef76SBarry Smith } 62543a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 62638597bd4SSatish Balay ierr = (*a->A->ops.mult)(a->A,xx,yy); CHKERRQ(ierr); 62743a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 62838597bd4SSatish Balay ierr = (*a->B->ops.multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr); 6291eb62cbbSBarry Smith return 0; 6301eb62cbbSBarry Smith } 6311eb62cbbSBarry Smith 6325615d1e5SSatish Balay #undef __FUNC__ 6335615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ" 6348f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 635da3a660dSBarry Smith { 636416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 637da3a660dSBarry Smith int ierr; 63843a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 63938597bd4SSatish Balay ierr = (*a->A->ops.multadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 64043a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 64138597bd4SSatish Balay ierr = (*a->B->ops.multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr); 642da3a660dSBarry Smith return 0; 643da3a660dSBarry Smith } 644da3a660dSBarry Smith 6455615d1e5SSatish Balay #undef __FUNC__ 6465615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ" 6478f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy) 648da3a660dSBarry Smith { 649416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 650da3a660dSBarry Smith int ierr; 651da3a660dSBarry Smith 652da3a660dSBarry Smith /* do nondiagonal part */ 65338597bd4SSatish Balay ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 654da3a660dSBarry Smith /* send it on its way */ 655537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 656da3a660dSBarry Smith /* do local part */ 65738597bd4SSatish Balay ierr = (*a->A->ops.multtrans)(a->A,xx,yy); CHKERRQ(ierr); 658da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 659da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 660da3a660dSBarry Smith /* but is not perhaps always true. */ 661537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 662da3a660dSBarry Smith return 0; 663da3a660dSBarry Smith } 664da3a660dSBarry Smith 6655615d1e5SSatish Balay #undef __FUNC__ 6665615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ" 6678f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 668da3a660dSBarry Smith { 669416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 670da3a660dSBarry Smith int ierr; 671da3a660dSBarry Smith 672da3a660dSBarry Smith /* do nondiagonal part */ 67338597bd4SSatish Balay ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 674da3a660dSBarry Smith /* send it on its way */ 675537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 676da3a660dSBarry Smith /* do local part */ 67738597bd4SSatish Balay ierr = (*a->A->ops.multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 678da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 679da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 680da3a660dSBarry Smith /* but is not perhaps always true. */ 6810a198c4cSBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 682da3a660dSBarry Smith return 0; 683da3a660dSBarry Smith } 684da3a660dSBarry Smith 6851eb62cbbSBarry Smith /* 6861eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 6871eb62cbbSBarry Smith diagonal block 6881eb62cbbSBarry Smith */ 6895615d1e5SSatish Balay #undef __FUNC__ 6905615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ" 6918f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v) 6921eb62cbbSBarry Smith { 693416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 69444cd7ae7SLois Curfman McInnes if (a->M != a->N) 695e3372554SBarry Smith SETERRQ(1,0,"Supports only square matrix where A->A is diag block"); 6965baf8537SBarry Smith if (a->rstart != a->cstart || a->rend != a->cend) { 697e3372554SBarry Smith SETERRQ(1,0,"row partition must equal col partition"); } 698416022c9SBarry Smith return MatGetDiagonal(a->A,v); 6991eb62cbbSBarry Smith } 7001eb62cbbSBarry Smith 7015615d1e5SSatish Balay #undef __FUNC__ 7025615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ" 7038f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A) 704052efed2SBarry Smith { 705052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 706052efed2SBarry Smith int ierr; 707052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 708052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 709052efed2SBarry Smith return 0; 710052efed2SBarry Smith } 711052efed2SBarry Smith 7125615d1e5SSatish Balay #undef __FUNC__ 7135eea60f9SBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" /* ADIC Ignore */ 7148f6be9afSLois Curfman McInnes int MatDestroy_MPIAIJ(PetscObject obj) 7151eb62cbbSBarry Smith { 7161eb62cbbSBarry Smith Mat mat = (Mat) obj; 71744a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 7181eb62cbbSBarry Smith int ierr; 71983e2fdc7SBarry Smith 720a5a9c739SBarry Smith #if defined(PETSC_LOG) 7216e35fa57SLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N); 722a5a9c739SBarry Smith #endif 72383e2fdc7SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 7240452661fSBarry Smith PetscFree(aij->rowners); 72578b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 72678b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 7270452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 7280452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 7291eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 730a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 7317a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 7320452661fSBarry Smith PetscFree(aij); 73390f02eecSBarry Smith if (mat->mapping) { 73490f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 73590f02eecSBarry Smith } 736a5a9c739SBarry Smith PLogObjectDestroy(mat); 7370452661fSBarry Smith PetscHeaderDestroy(mat); 7381eb62cbbSBarry Smith return 0; 7391eb62cbbSBarry Smith } 740ee50ffe9SBarry Smith 7415615d1e5SSatish Balay #undef __FUNC__ 7425eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" /* ADIC Ignore */ 7438f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 7441eb62cbbSBarry Smith { 745416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 746416022c9SBarry Smith int ierr; 747416022c9SBarry Smith 74817699dbbSLois Curfman McInnes if (aij->size == 1) { 749416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 750416022c9SBarry Smith } 751e3372554SBarry Smith else SETERRQ(1,0,"Only uniprocessor output supported"); 752416022c9SBarry Smith return 0; 753416022c9SBarry Smith } 754416022c9SBarry Smith 7555615d1e5SSatish Balay #undef __FUNC__ 7565eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab" /* ADIC Ignore */ 7578f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer) 758416022c9SBarry Smith { 75944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 760dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 761a56f8943SBarry Smith int ierr, format,shift = C->indexshift,rank; 762d636dbe3SBarry Smith FILE *fd; 76319bcc07fSBarry Smith ViewerType vtype; 764416022c9SBarry Smith 76519bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 76619bcc07fSBarry Smith if (vtype == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) { 76790ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 7680a198c4cSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 7694e220ebcSLois Curfman McInnes MatInfo info; 7704e220ebcSLois Curfman McInnes int flg; 771a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 77290ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 7734e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 77495e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 77577c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 77695e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 7774e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 77895e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 7794e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 7804e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->A,MAT_LOCAL,&info); 7814e220ebcSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used); 7824e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->B,MAT_LOCAL,&info); 7834e220ebcSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used); 784a56f8943SBarry Smith fflush(fd); 78577c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 786a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 787416022c9SBarry Smith return 0; 788416022c9SBarry Smith } 7890a198c4cSBarry Smith else if (format == VIEWER_FORMAT_ASCII_INFO) { 79008480c60SBarry Smith return 0; 79108480c60SBarry Smith } 792416022c9SBarry Smith } 793416022c9SBarry Smith 79419bcc07fSBarry Smith if (vtype == DRAW_VIEWER) { 79519bcc07fSBarry Smith Draw draw; 79619bcc07fSBarry Smith PetscTruth isnull; 79719bcc07fSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 79819bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 79919bcc07fSBarry Smith } 80019bcc07fSBarry Smith 80119bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 80290ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 80377c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 804d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 80517699dbbSLois Curfman McInnes aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 8061eb62cbbSBarry Smith aij->cend); 80778b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 80878b31e54SBarry Smith ierr = MatView(aij->B,viewer); CHKERRQ(ierr); 809d13ab20cSBarry Smith fflush(fd); 81077c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 811d13ab20cSBarry Smith } 812416022c9SBarry Smith else { 813a56f8943SBarry Smith int size = aij->size; 814a56f8943SBarry Smith rank = aij->rank; 81517699dbbSLois Curfman McInnes if (size == 1) { 81678b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 81795373324SBarry Smith } 81895373324SBarry Smith else { 81995373324SBarry Smith /* assemble the entire matrix onto first processor. */ 82095373324SBarry Smith Mat A; 821ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 8222eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 82395373324SBarry Smith Scalar *a; 8242ee70a88SLois Curfman McInnes 82517699dbbSLois Curfman McInnes if (!rank) { 826b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 827c451ab25SLois Curfman McInnes CHKERRQ(ierr); 82895373324SBarry Smith } 82995373324SBarry Smith else { 830b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 831c451ab25SLois Curfman McInnes CHKERRQ(ierr); 83295373324SBarry Smith } 833464493b3SBarry Smith PLogObjectParent(mat,A); 834416022c9SBarry Smith 83595373324SBarry Smith /* copy over the A part */ 836ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 8372ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 83895373324SBarry Smith row = aij->rstart; 839dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 84095373324SBarry Smith for ( i=0; i<m; i++ ) { 841416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 84295373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 84395373324SBarry Smith } 8442ee70a88SLois Curfman McInnes aj = Aloc->j; 845dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 84695373324SBarry Smith 84795373324SBarry Smith /* copy over the B part */ 848ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 8492ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 85095373324SBarry Smith row = aij->rstart; 8510452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 852dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 85395373324SBarry Smith for ( i=0; i<m; i++ ) { 854416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 85595373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 85695373324SBarry Smith } 8570452661fSBarry Smith PetscFree(ct); 8586d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8596d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 86017699dbbSLois Curfman McInnes if (!rank) { 86178b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 86295373324SBarry Smith } 86378b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 86495373324SBarry Smith } 86595373324SBarry Smith } 8661eb62cbbSBarry Smith return 0; 8671eb62cbbSBarry Smith } 8681eb62cbbSBarry Smith 8695615d1e5SSatish Balay #undef __FUNC__ 8705eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ" /* ADIC Ignore */ 8718f6be9afSLois Curfman McInnes int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 872416022c9SBarry Smith { 873416022c9SBarry Smith Mat mat = (Mat) obj; 874416022c9SBarry Smith int ierr; 87519bcc07fSBarry Smith ViewerType vtype; 876416022c9SBarry Smith 87719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 87819bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER || 87919bcc07fSBarry Smith vtype == DRAW_VIEWER || vtype == MATLAB_VIEWER) { 880d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 881416022c9SBarry Smith } 88219bcc07fSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 883416022c9SBarry Smith return MatView_MPIAIJ_Binary(mat,viewer); 884416022c9SBarry Smith } 885416022c9SBarry Smith return 0; 886416022c9SBarry Smith } 887416022c9SBarry Smith 8881eb62cbbSBarry Smith /* 8891eb62cbbSBarry Smith This has to provide several versions. 8901eb62cbbSBarry Smith 8911eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 8921eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 893d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 8941eb62cbbSBarry Smith */ 8955615d1e5SSatish Balay #undef __FUNC__ 8965615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ" 8978f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 898dbb450caSBarry Smith double fshift,int its,Vec xx) 8998a729477SBarry Smith { 90044a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 901d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 902ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 903c16cb8f2SBarry Smith Scalar *b,*x,*xs,*ls,d,*v,sum; 9046abc6512SBarry Smith int ierr,*idx, *diag; 905416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 9068a729477SBarry Smith 907d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 908dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 909dbb450caSBarry Smith ls = ls + shift; 91083e2fdc7SBarry Smith if (!A->diag) {ierr = MatMarkDiag_SeqAIJ(AA); CHKERRQ(ierr);} 911d6dfbf8fSBarry Smith diag = A->diag; 912c16cb8f2SBarry Smith if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 913da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 91438597bd4SSatish Balay return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx); 915da3a660dSBarry Smith } 91643a90d84SBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 91778b31e54SBarry Smith CHKERRQ(ierr); 91843a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 91978b31e54SBarry Smith CHKERRQ(ierr); 920d6dfbf8fSBarry Smith while (its--) { 921d6dfbf8fSBarry Smith /* go down through the rows */ 922d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 923d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 924dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 925dbb450caSBarry Smith v = A->a + A->i[i] + shift; 926d6dfbf8fSBarry Smith sum = b[i]; 927d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 928dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 929d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 930dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 931dbb450caSBarry Smith v = B->a + B->i[i] + shift; 932d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 93355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 934d6dfbf8fSBarry Smith } 935d6dfbf8fSBarry Smith /* come up through the rows */ 936d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 937d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 938dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 939dbb450caSBarry Smith v = A->a + A->i[i] + shift; 940d6dfbf8fSBarry Smith sum = b[i]; 941d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 942dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 943d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 944dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 945dbb450caSBarry Smith v = B->a + B->i[i] + shift; 946d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 94755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 948d6dfbf8fSBarry Smith } 949d6dfbf8fSBarry Smith } 950d6dfbf8fSBarry Smith } 951d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 952da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 95338597bd4SSatish Balay return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx); 954da3a660dSBarry Smith } 95543a90d84SBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 95678b31e54SBarry Smith CHKERRQ(ierr); 95743a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 95878b31e54SBarry Smith CHKERRQ(ierr); 959d6dfbf8fSBarry Smith while (its--) { 960d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 961d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 962dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 963dbb450caSBarry Smith v = A->a + A->i[i] + shift; 964d6dfbf8fSBarry Smith sum = b[i]; 965d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 966dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 967d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 968dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 969dbb450caSBarry Smith v = B->a + B->i[i] + shift; 970d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 97155a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 972d6dfbf8fSBarry Smith } 973d6dfbf8fSBarry Smith } 974d6dfbf8fSBarry Smith } 975d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 976da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 97738597bd4SSatish Balay return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx); 978da3a660dSBarry Smith } 97943a90d84SBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 98078b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 98143a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 98278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 983d6dfbf8fSBarry Smith while (its--) { 984d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 985d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 986dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 987dbb450caSBarry Smith v = A->a + A->i[i] + shift; 988d6dfbf8fSBarry Smith sum = b[i]; 989d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 990dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 991d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 992dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 993dbb450caSBarry Smith v = B->a + B->i[i] + shift; 994d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 99555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 996d6dfbf8fSBarry Smith } 997d6dfbf8fSBarry Smith } 998d6dfbf8fSBarry Smith } 999c16cb8f2SBarry Smith else { 1000e3372554SBarry Smith SETERRQ(1,0,"Option not supported"); 1001c16cb8f2SBarry Smith } 10028a729477SBarry Smith return 0; 10038a729477SBarry Smith } 1004a66be287SLois Curfman McInnes 10055615d1e5SSatish Balay #undef __FUNC__ 10065eea60f9SBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" /* ADIC Ignore */ 10078f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info) 1008a66be287SLois Curfman McInnes { 1009a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1010a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 10114e220ebcSLois Curfman McInnes int ierr; 10124e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 1013a66be287SLois Curfman McInnes 10144e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 10154e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 10164e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 10174e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 10184e220ebcSLois Curfman McInnes info->block_size = 1.0; 10194e220ebcSLois Curfman McInnes ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr); 10204e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 10214e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 10224e220ebcSLois Curfman McInnes ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr); 10234e220ebcSLois Curfman McInnes isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded; 10244e220ebcSLois Curfman McInnes isend[3] += info->memory; isend[4] += info->mallocs; 1025a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 10264e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 10274e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 10284e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 10294e220ebcSLois Curfman McInnes info->memory = isend[3]; 10304e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 1031a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 10324e220ebcSLois Curfman McInnes MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm); 10334e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10344e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10354e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10364e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10374e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1038a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 10394e220ebcSLois Curfman McInnes MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm); 10404e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10414e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10424e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10434e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10444e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1045a66be287SLois Curfman McInnes } 10464e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 10474e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10484e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10494e220ebcSLois Curfman McInnes 1050a66be287SLois Curfman McInnes return 0; 1051a66be287SLois Curfman McInnes } 1052a66be287SLois Curfman McInnes 10535615d1e5SSatish Balay #undef __FUNC__ 10545eea60f9SBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" /* ADIC Ignore */ 10558f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op) 1056c74985f6SBarry Smith { 1057c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1058c74985f6SBarry Smith 10596d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 10606d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 10616da5968aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 1062c2653b3dSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 106396854ed6SLois Curfman McInnes op == MAT_NEW_NONZERO_ALLOCATION_ERROR || 1064c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR) { 1065b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 1066b1fbbac0SLois Curfman McInnes MatSetOption(a->B,op); 1067b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 1068aeafbbfcSLois Curfman McInnes a->roworiented = 1; 1069c0bbcb79SLois Curfman McInnes MatSetOption(a->A,op); 1070c0bbcb79SLois Curfman McInnes MatSetOption(a->B,op); 1071b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 10726da5968aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 10736d4a8577SBarry Smith op == MAT_SYMMETRIC || 10746d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 10756d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 107694a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_MPIAIJ:Option ignored\n"); 10776d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) { 10784b0e389bSBarry Smith a->roworiented = 0; 10794b0e389bSBarry Smith MatSetOption(a->A,op); 10804b0e389bSBarry Smith MatSetOption(a->B,op); 10812b362799SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 108290f02eecSBarry Smith a->donotstash = 1; 108390f02eecSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) 1084e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");} 1085c0bbcb79SLois Curfman McInnes else 1086e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"unknown option");} 1087c74985f6SBarry Smith return 0; 1088c74985f6SBarry Smith } 1089c74985f6SBarry Smith 10905615d1e5SSatish Balay #undef __FUNC__ 10915eea60f9SBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" /* ADIC Ignore */ 10928f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1093c74985f6SBarry Smith { 109444a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1095c74985f6SBarry Smith *m = mat->M; *n = mat->N; 1096c74985f6SBarry Smith return 0; 1097c74985f6SBarry Smith } 1098c74985f6SBarry Smith 10995615d1e5SSatish Balay #undef __FUNC__ 11005eea60f9SBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" /* ADIC Ignore */ 11018f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1102c74985f6SBarry Smith { 110344a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1104b7c46309SBarry Smith *m = mat->m; *n = mat->N; 1105c74985f6SBarry Smith return 0; 1106c74985f6SBarry Smith } 1107c74985f6SBarry Smith 11085615d1e5SSatish Balay #undef __FUNC__ 11095eea60f9SBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" /* ADIC Ignore */ 11108f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1111c74985f6SBarry Smith { 111244a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1113c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 1114c74985f6SBarry Smith return 0; 1115c74985f6SBarry Smith } 1116c74985f6SBarry Smith 11176d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11186d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11196d84be18SBarry Smith 11205615d1e5SSatish Balay #undef __FUNC__ 11215615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ" 11226d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 112339e00950SLois Curfman McInnes { 1124154123eaSLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 112570f0671dSBarry Smith Scalar *vworkA, *vworkB, **pvA, **pvB,*v_p; 1126154123eaSLois Curfman McInnes int i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart; 1127154123eaSLois Curfman McInnes int nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend; 112870f0671dSBarry Smith int *cmap, *idx_p; 112939e00950SLois Curfman McInnes 1130e3372554SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,0,"Already active"); 11317a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11327a0afa10SBarry Smith 113370f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11347a0afa10SBarry Smith /* 11357a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11367a0afa10SBarry Smith */ 11377a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 1138c16cb8f2SBarry Smith int max = 1,m = mat->m,tmp; 1139c16cb8f2SBarry Smith for ( i=0; i<m; i++ ) { 11407a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11417a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11427a0afa10SBarry Smith } 11437a0afa10SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar))); 11447a0afa10SBarry Smith CHKPTRQ(mat->rowvalues); 11457a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11467a0afa10SBarry Smith } 11477a0afa10SBarry Smith 1148e3372554SBarry Smith if (row < rstart || row >= rend) SETERRQ(1,0,"Only local rows") 1149abc0e9e4SLois Curfman McInnes lrow = row - rstart; 115039e00950SLois Curfman McInnes 1151154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1152154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1153154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 115438597bd4SSatish Balay ierr = (*mat->A->ops.getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 115538597bd4SSatish Balay ierr = (*mat->B->ops.getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1156154123eaSLois Curfman McInnes nztot = nzA + nzB; 1157154123eaSLois Curfman McInnes 115870f0671dSBarry Smith cmap = mat->garray; 1159154123eaSLois Curfman McInnes if (v || idx) { 1160154123eaSLois Curfman McInnes if (nztot) { 1161154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 116270f0671dSBarry Smith int imark = -1; 1163154123eaSLois Curfman McInnes if (v) { 116470f0671dSBarry Smith *v = v_p = mat->rowvalues; 116539e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 116670f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1167154123eaSLois Curfman McInnes else break; 1168154123eaSLois Curfman McInnes } 1169154123eaSLois Curfman McInnes imark = i; 117070f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 117170f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1172154123eaSLois Curfman McInnes } 1173154123eaSLois Curfman McInnes if (idx) { 117470f0671dSBarry Smith *idx = idx_p = mat->rowindices; 117570f0671dSBarry Smith if (imark > -1) { 117670f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 117770f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 117870f0671dSBarry Smith } 117970f0671dSBarry Smith } else { 1180154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 118170f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1182154123eaSLois Curfman McInnes else break; 1183154123eaSLois Curfman McInnes } 1184154123eaSLois Curfman McInnes imark = i; 118570f0671dSBarry Smith } 118670f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 118770f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 118839e00950SLois Curfman McInnes } 118939e00950SLois Curfman McInnes } 11901ca473b0SSatish Balay else { 11911ca473b0SSatish Balay if (idx) *idx = 0; 11921ca473b0SSatish Balay if (v) *v = 0; 11931ca473b0SSatish Balay } 1194154123eaSLois Curfman McInnes } 119539e00950SLois Curfman McInnes *nz = nztot; 119638597bd4SSatish Balay ierr = (*mat->A->ops.restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 119738597bd4SSatish Balay ierr = (*mat->B->ops.restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 119839e00950SLois Curfman McInnes return 0; 119939e00950SLois Curfman McInnes } 120039e00950SLois Curfman McInnes 12015615d1e5SSatish Balay #undef __FUNC__ 12025eea60f9SBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" /* ADIC Ignore */ 12036d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 120439e00950SLois Curfman McInnes { 12057a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 12067a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 1207e3372554SBarry Smith SETERRQ(1,0,"MatGetRow not called"); 12087a0afa10SBarry Smith } 12097a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 121039e00950SLois Curfman McInnes return 0; 121139e00950SLois Curfman McInnes } 121239e00950SLois Curfman McInnes 12135615d1e5SSatish Balay #undef __FUNC__ 12145615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ" 12158f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1216855ac2c5SLois Curfman McInnes { 1217855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1218ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1219416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1220416022c9SBarry Smith double sum = 0.0; 122104ca555eSLois Curfman McInnes Scalar *v; 122204ca555eSLois Curfman McInnes 122317699dbbSLois Curfman McInnes if (aij->size == 1) { 122414183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 122537fa93a5SLois Curfman McInnes } else { 122604ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 122704ca555eSLois Curfman McInnes v = amat->a; 122804ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 122904ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 123004ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 123104ca555eSLois Curfman McInnes #else 123204ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 123304ca555eSLois Curfman McInnes #endif 123404ca555eSLois Curfman McInnes } 123504ca555eSLois Curfman McInnes v = bmat->a; 123604ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 123704ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 123804ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 123904ca555eSLois Curfman McInnes #else 124004ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 124104ca555eSLois Curfman McInnes #endif 124204ca555eSLois Curfman McInnes } 12436d84be18SBarry Smith MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm); 124404ca555eSLois Curfman McInnes *norm = sqrt(*norm); 124504ca555eSLois Curfman McInnes } 124604ca555eSLois Curfman McInnes else if (type == NORM_1) { /* max column norm */ 124704ca555eSLois Curfman McInnes double *tmp, *tmp2; 124804ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 1249758f045eSSatish Balay tmp = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp); 1250758f045eSSatish Balay tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2); 1251cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 125204ca555eSLois Curfman McInnes *norm = 0.0; 125304ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 125404ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1255579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 125604ca555eSLois Curfman McInnes } 125704ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 125804ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1259579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 126004ca555eSLois Curfman McInnes } 12616d84be18SBarry Smith MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm); 126204ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 126304ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 126404ca555eSLois Curfman McInnes } 12650452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 126604ca555eSLois Curfman McInnes } 126704ca555eSLois Curfman McInnes else if (type == NORM_INFINITY) { /* max row norm */ 1268515d9167SLois Curfman McInnes double ntemp = 0.0; 126904ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1270dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 127104ca555eSLois Curfman McInnes sum = 0.0; 127204ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1273cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 127404ca555eSLois Curfman McInnes } 1275dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 127604ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1277cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 127804ca555eSLois Curfman McInnes } 1279515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 128004ca555eSLois Curfman McInnes } 12816d84be18SBarry Smith MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm); 128204ca555eSLois Curfman McInnes } 128304ca555eSLois Curfman McInnes else { 1284e3372554SBarry Smith SETERRQ(1,0,"No support for two norm"); 128504ca555eSLois Curfman McInnes } 128637fa93a5SLois Curfman McInnes } 1287855ac2c5SLois Curfman McInnes return 0; 1288855ac2c5SLois Curfman McInnes } 1289855ac2c5SLois Curfman McInnes 12905615d1e5SSatish Balay #undef __FUNC__ 12915615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ" 12928f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1293b7c46309SBarry Smith { 1294b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1295dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1296416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1297416022c9SBarry Smith Mat B; 1298b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 1299b7c46309SBarry Smith Scalar *array; 1300b7c46309SBarry Smith 13013638b69dSLois Curfman McInnes if (matout == PETSC_NULL && M != N) 1302e3372554SBarry Smith SETERRQ(1,0,"Square matrix only for in-place"); 1303b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,PETSC_NULL,0, 1304b4fd4287SBarry Smith PETSC_NULL,&B); CHKERRQ(ierr); 1305b7c46309SBarry Smith 1306b7c46309SBarry Smith /* copy over the A part */ 1307ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1308b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1309b7c46309SBarry Smith row = a->rstart; 1310dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1311b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1312416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1313b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1314b7c46309SBarry Smith } 1315b7c46309SBarry Smith aj = Aloc->j; 13164af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1317b7c46309SBarry Smith 1318b7c46309SBarry Smith /* copy over the B part */ 1319ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1320b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1321b7c46309SBarry Smith row = a->rstart; 13220452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1323dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1324b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1325416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1326b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1327b7c46309SBarry Smith } 13280452661fSBarry Smith PetscFree(ct); 13296d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13306d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13313638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13320de55854SLois Curfman McInnes *matout = B; 13330de55854SLois Curfman McInnes } else { 13340de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13350452661fSBarry Smith PetscFree(a->rowners); 13360de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13370de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 13380452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 13390452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13400de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1341a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13420452661fSBarry Smith PetscFree(a); 1343f09e8eb9SSatish Balay PetscMemcpy(A,B,sizeof(struct _p_Mat)); 13440452661fSBarry Smith PetscHeaderDestroy(B); 13450de55854SLois Curfman McInnes } 1346b7c46309SBarry Smith return 0; 1347b7c46309SBarry Smith } 1348b7c46309SBarry Smith 13495615d1e5SSatish Balay #undef __FUNC__ 13505615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ" 13514b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr) 1352a008b906SSatish Balay { 13534b967eb1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 13544b967eb1SSatish Balay Mat a = aij->A, b = aij->B; 1355a008b906SSatish Balay int ierr,s1,s2,s3; 1356a008b906SSatish Balay 13574b967eb1SSatish Balay ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr); 13584b967eb1SSatish Balay if (rr) { 13594b967eb1SSatish Balay s3 = aij->n; 13604b967eb1SSatish Balay VecGetLocalSize_Fast(rr,s1); 1361e3372554SBarry Smith if (s1!=s3) SETERRQ(1,0,"right vector non-conforming local size"); 13624b967eb1SSatish Balay /* Overlap communication with computation. */ 136343a90d84SBarry Smith ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1364a008b906SSatish Balay } 13654b967eb1SSatish Balay if (ll) { 13664b967eb1SSatish Balay VecGetLocalSize_Fast(ll,s1); 1367e3372554SBarry Smith if (s1!=s2) SETERRQ(1,0,"left vector non-conforming local size"); 1368c351683dSSatish Balay ierr = (*b->ops.diagonalscale)(b,ll,0); CHKERRQ(ierr); 13694b967eb1SSatish Balay } 13704b967eb1SSatish Balay /* scale the diagonal block */ 1371c351683dSSatish Balay ierr = (*a->ops.diagonalscale)(a,ll,rr); CHKERRQ(ierr); 13724b967eb1SSatish Balay 13734b967eb1SSatish Balay if (rr) { 13744b967eb1SSatish Balay /* Do a scatter end and then right scale the off-diagonal block */ 137543a90d84SBarry Smith ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1376c351683dSSatish Balay ierr = (*b->ops.diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr); 13774b967eb1SSatish Balay } 13784b967eb1SSatish Balay 1379a008b906SSatish Balay return 0; 1380a008b906SSatish Balay } 1381a008b906SSatish Balay 1382a008b906SSatish Balay 1383682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat); 13845615d1e5SSatish Balay #undef __FUNC__ 13855eea60f9SBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" /* ADIC Ignore */ 13868f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A) 1387682d7d0cSBarry Smith { 1388682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1389682d7d0cSBarry Smith 1390682d7d0cSBarry Smith if (!a->rank) return MatPrintHelp_SeqAIJ(a->A); 1391682d7d0cSBarry Smith else return 0; 1392682d7d0cSBarry Smith } 1393682d7d0cSBarry Smith 13945615d1e5SSatish Balay #undef __FUNC__ 13955eea60f9SBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" /* ADIC Ignore */ 13968f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs) 13975a838052SSatish Balay { 13985a838052SSatish Balay *bs = 1; 13995a838052SSatish Balay return 0; 14005a838052SSatish Balay } 14015615d1e5SSatish Balay #undef __FUNC__ 14025eea60f9SBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" /* ADIC Ignore */ 14038f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A) 1404bb5a7306SBarry Smith { 1405bb5a7306SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1406bb5a7306SBarry Smith int ierr; 1407bb5a7306SBarry Smith ierr = MatSetUnfactored(a->A); CHKERRQ(ierr); 1408bb5a7306SBarry Smith return 0; 1409bb5a7306SBarry Smith } 1410bb5a7306SBarry Smith 14118f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int); 14122f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 14130a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 14140a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **); 14158a729477SBarry Smith /* -------------------------------------------------------------------*/ 14162ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ, 141739e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 141844a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 141944a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 142036ce4990SBarry Smith 0,0, 142136ce4990SBarry Smith 0,0, 142236ce4990SBarry Smith 0,0, 142344a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1424b7c46309SBarry Smith MatTranspose_MPIAIJ, 1425a66be287SLois Curfman McInnes MatGetInfo_MPIAIJ,0, 1426a008b906SSatish Balay MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ, 1427ee50ffe9SBarry Smith MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ, 14281eb62cbbSBarry Smith 0, 1429299609e3SLois Curfman McInnes MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ, 143036ce4990SBarry Smith 0,0,0,0, 1431d1710a03SLois Curfman McInnes MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ, 143236ce4990SBarry Smith 0,0, 143394a9d846SBarry Smith 0,0,MatConvertSameType_MPIAIJ,0,0, 1434b49de8d1SLois Curfman McInnes 0,0,0, 1435598137ffSSatish Balay MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0, 1436052efed2SBarry Smith MatPrintHelp_MPIAIJ, 14373b2fbd54SBarry Smith MatScale_MPIAIJ,0,0,0, 14380a198c4cSBarry Smith MatGetBlockSize_MPIAIJ,0,0,0,0, 1439bb5a7306SBarry Smith MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ}; 144036ce4990SBarry Smith 14418a729477SBarry Smith 14425615d1e5SSatish Balay #undef __FUNC__ 14435615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ" 14441987afe7SBarry Smith /*@C 1445ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 14463a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 14473a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 14483a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 14493a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 14508a729477SBarry Smith 14518a729477SBarry Smith Input Parameters: 14521eb62cbbSBarry Smith . comm - MPI communicator 14537d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 145492e8d321SLois Curfman McInnes This value should be the same as the local size used in creating the 145592e8d321SLois Curfman McInnes y vector for the matrix-vector product y = Ax. 14561a3896d6SBarry Smith . n - This value should be the same as the local size used in creating the 14571a3896d6SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 14581a3896d6SBarry Smith calculated if N is given) For square matrices n is almost always m. 14597d3e4905SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 146092e8d321SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1461ab693e5aSLois Curfman McInnes . d_nz - number of nonzeros per row in diagonal portion of local submatrix 1462ff756334SLois Curfman McInnes (same for all local rows) 14632bd5e0b2SLois Curfman McInnes . d_nzz - array containing the number of nonzeros in the various rows of the 146492e8d321SLois Curfman McInnes diagonal portion of the local submatrix (possibly different for each row) 14652bd5e0b2SLois Curfman McInnes or PETSC_NULL. You must leave room for the diagonal entry even if 14662bd5e0b2SLois Curfman McInnes it is zero. 14672bd5e0b2SLois Curfman McInnes . o_nz - number of nonzeros per row in the off-diagonal portion of local 1468ab693e5aSLois Curfman McInnes submatrix (same for all local rows). 14692bd5e0b2SLois Curfman McInnes . o_nzz - array containing the number of nonzeros in the various rows of the 14702bd5e0b2SLois Curfman McInnes off-diagonal portion of the local submatrix (possibly different for 14712bd5e0b2SLois Curfman McInnes each row) or PETSC_NULL. 14728a729477SBarry Smith 1473ff756334SLois Curfman McInnes Output Parameter: 147444cd7ae7SLois Curfman McInnes . A - the matrix 14758a729477SBarry Smith 1476ff756334SLois Curfman McInnes Notes: 1477ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1478ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 14790002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 14800002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1481ff756334SLois Curfman McInnes 1482ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1483ff756334SLois Curfman McInnes (possibly both). 1484ff756334SLois Curfman McInnes 14855511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 14865511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 14875511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 14885511cfe3SLois Curfman McInnes 14895511cfe3SLois Curfman McInnes Options Database Keys: 14905511cfe3SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 14916e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 14926e62573dSLois Curfman McInnes $ (max limit=5) 14936e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 14946e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 14956e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 14965511cfe3SLois Curfman McInnes 1497e0245417SLois Curfman McInnes Storage Information: 1498e0245417SLois Curfman McInnes For a square global matrix we define each processor's diagonal portion 1499e0245417SLois Curfman McInnes to be its local rows and the corresponding columns (a square submatrix); 1500e0245417SLois Curfman McInnes each processor's off-diagonal portion encompasses the remainder of the 1501e0245417SLois Curfman McInnes local matrix (a rectangular submatrix). 1502e0245417SLois Curfman McInnes 1503e0245417SLois Curfman McInnes The user can specify preallocated storage for the diagonal part of 15045ace5be8SLois Curfman McInnes the local submatrix with either d_nz or d_nnz (not both). Set 15055ace5be8SLois Curfman McInnes d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic 15065ace5be8SLois Curfman McInnes memory allocation. Likewise, specify preallocated storage for the 15075ace5be8SLois Curfman McInnes off-diagonal part of the local submatrix with o_nz or o_nnz (not both). 1508ff756334SLois Curfman McInnes 15095511cfe3SLois Curfman McInnes Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In 15105511cfe3SLois Curfman McInnes the figure below we depict these three local rows and all columns (0-11). 15112191d07cSBarry Smith 1512b810aeb4SBarry Smith $ 0 1 2 3 4 5 6 7 8 9 10 11 1513b810aeb4SBarry Smith $ ------------------- 15145511cfe3SLois Curfman McInnes $ row 3 | o o o d d d o o o o o o 15155511cfe3SLois Curfman McInnes $ row 4 | o o o d d d o o o o o o 15165511cfe3SLois Curfman McInnes $ row 5 | o o o d d d o o o o o o 15175511cfe3SLois Curfman McInnes $ ------------------- 1518b810aeb4SBarry Smith $ 1519b810aeb4SBarry Smith 15205511cfe3SLois Curfman McInnes Thus, any entries in the d locations are stored in the d (diagonal) 15215511cfe3SLois Curfman McInnes submatrix, and any entries in the o locations are stored in the 15225511cfe3SLois Curfman McInnes o (off-diagonal) submatrix. Note that the d and the o submatrices are 15235511cfe3SLois Curfman McInnes stored simply in the MATSEQAIJ format for compressed row storage. 15245511cfe3SLois Curfman McInnes 15255511cfe3SLois Curfman McInnes Now d_nz should indicate the number of nonzeros per row in the d matrix, 15265511cfe3SLois Curfman McInnes and o_nz should indicate the number of nonzeros per row in the o matrix. 15275511cfe3SLois Curfman McInnes In general, for PDE problems in which most nonzeros are near the diagonal, 15283d323bbdSBarry Smith one expects d_nz >> o_nz. For large problems you MUST preallocate memory 152992e8d321SLois Curfman McInnes or you will get TERRIBLE performance; see the users' manual chapter on 15306da5968aSLois Curfman McInnes matrices. 15313a511b96SLois Curfman McInnes 1532dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1533ff756334SLois Curfman McInnes 1534fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 15358a729477SBarry Smith @*/ 15361eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 153744cd7ae7SLois Curfman McInnes int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A) 15388a729477SBarry Smith { 153944cd7ae7SLois Curfman McInnes Mat B; 154044cd7ae7SLois Curfman McInnes Mat_MPIAIJ *b; 154136ce4990SBarry Smith int ierr, i,sum[2],work[2],size; 1542416022c9SBarry Smith 15433914022bSBarry Smith MPI_Comm_size(comm,&size); 15443914022bSBarry Smith if (size == 1) { 15453914022bSBarry Smith if (M == PETSC_DECIDE) M = m; 15463914022bSBarry Smith if (N == PETSC_DECIDE) N = n; 15473914022bSBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,d_nz,d_nnz,A); CHKERRQ(ierr); 15483914022bSBarry Smith return 0; 15493914022bSBarry Smith } 15503914022bSBarry Smith 155144cd7ae7SLois Curfman McInnes *A = 0; 1552f09e8eb9SSatish Balay PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATMPIAIJ,comm); 155344cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 155444cd7ae7SLois Curfman McInnes B->data = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b); 155544cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_MPIAIJ)); 155644cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 155744cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_MPIAIJ; 155844cd7ae7SLois Curfman McInnes B->view = MatView_MPIAIJ; 155944cd7ae7SLois Curfman McInnes B->factor = 0; 156044cd7ae7SLois Curfman McInnes B->assembled = PETSC_FALSE; 156190f02eecSBarry Smith B->mapping = 0; 1562d6dfbf8fSBarry Smith 156347794344SBarry Smith B->insertmode = NOT_SET_VALUES; 1564*9eb4d147SSatish Balay b->size = size; 156544cd7ae7SLois Curfman McInnes MPI_Comm_rank(comm,&b->rank); 15661eb62cbbSBarry Smith 1567b4fd4287SBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) 1568e3372554SBarry Smith SETERRQ(1,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz"); 15691987afe7SBarry Smith 1570dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE || N == PETSC_DECIDE) { 15711eb62cbbSBarry Smith work[0] = m; work[1] = n; 1572d65a2f8fSBarry Smith MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm ); 1573dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE) M = sum[0]; 1574dbd7a890SLois Curfman McInnes if (N == PETSC_DECIDE) N = sum[1]; 15751eb62cbbSBarry Smith } 157644cd7ae7SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);} 157744cd7ae7SLois Curfman McInnes if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);} 157844cd7ae7SLois Curfman McInnes b->m = m; B->m = m; 157944cd7ae7SLois Curfman McInnes b->n = n; B->n = n; 158044cd7ae7SLois Curfman McInnes b->N = N; B->N = N; 158144cd7ae7SLois Curfman McInnes b->M = M; B->M = M; 15821eb62cbbSBarry Smith 15831eb62cbbSBarry Smith /* build local table of row and column ownerships */ 158444cd7ae7SLois Curfman McInnes b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners); 1585f09e8eb9SSatish Balay PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1586603f58a4SSatish Balay b->cowners = b->rowners + b->size + 2; 158744cd7ae7SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm); 158844cd7ae7SLois Curfman McInnes b->rowners[0] = 0; 158944cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 159044cd7ae7SLois Curfman McInnes b->rowners[i] += b->rowners[i-1]; 15918a729477SBarry Smith } 159244cd7ae7SLois Curfman McInnes b->rstart = b->rowners[b->rank]; 159344cd7ae7SLois Curfman McInnes b->rend = b->rowners[b->rank+1]; 159444cd7ae7SLois Curfman McInnes MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm); 159544cd7ae7SLois Curfman McInnes b->cowners[0] = 0; 159644cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 159744cd7ae7SLois Curfman McInnes b->cowners[i] += b->cowners[i-1]; 15988a729477SBarry Smith } 159944cd7ae7SLois Curfman McInnes b->cstart = b->cowners[b->rank]; 160044cd7ae7SLois Curfman McInnes b->cend = b->cowners[b->rank+1]; 16018a729477SBarry Smith 16025ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 1603029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr); 160444cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->A); 16057b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 1606029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr); 160744cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->B); 16088a729477SBarry Smith 16091eb62cbbSBarry Smith /* build cache for off array entries formed */ 161044cd7ae7SLois Curfman McInnes ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr); 161190f02eecSBarry Smith b->donotstash = 0; 161244cd7ae7SLois Curfman McInnes b->colmap = 0; 161344cd7ae7SLois Curfman McInnes b->garray = 0; 161444cd7ae7SLois Curfman McInnes b->roworiented = 1; 16158a729477SBarry Smith 16161eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 161744cd7ae7SLois Curfman McInnes b->lvec = 0; 161844cd7ae7SLois Curfman McInnes b->Mvctx = 0; 16198a729477SBarry Smith 16207a0afa10SBarry Smith /* stuff for MatGetRow() */ 162144cd7ae7SLois Curfman McInnes b->rowindices = 0; 162244cd7ae7SLois Curfman McInnes b->rowvalues = 0; 162344cd7ae7SLois Curfman McInnes b->getrowactive = PETSC_FALSE; 16247a0afa10SBarry Smith 162544cd7ae7SLois Curfman McInnes *A = B; 1626d6dfbf8fSBarry Smith return 0; 1627d6dfbf8fSBarry Smith } 1628c74985f6SBarry Smith 16295615d1e5SSatish Balay #undef __FUNC__ 16305615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ" 16318f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues) 1632d6dfbf8fSBarry Smith { 1633d6dfbf8fSBarry Smith Mat mat; 1634416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 1635a1b97e82SLois Curfman McInnes int ierr, len=0, flg; 1636d6dfbf8fSBarry Smith 1637416022c9SBarry Smith *newmat = 0; 1638f09e8eb9SSatish Balay PetscHeaderCreate(mat,_p_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm); 1639a5a9c739SBarry Smith PLogObjectCreate(mat); 16400452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1641416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 164244a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 164344a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1644d6dfbf8fSBarry Smith mat->factor = matin->factor; 1645c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1646d6dfbf8fSBarry Smith 164744cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 164844cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 164944cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 165044cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 1651d6dfbf8fSBarry Smith 1652416022c9SBarry Smith a->rstart = oldmat->rstart; 1653416022c9SBarry Smith a->rend = oldmat->rend; 1654416022c9SBarry Smith a->cstart = oldmat->cstart; 1655416022c9SBarry Smith a->cend = oldmat->cend; 165617699dbbSLois Curfman McInnes a->size = oldmat->size; 165717699dbbSLois Curfman McInnes a->rank = oldmat->rank; 165847794344SBarry Smith mat->insertmode = NOT_SET_VALUES; 1659bcd2baecSBarry Smith a->rowvalues = 0; 1660bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1661d6dfbf8fSBarry Smith 1662603f58a4SSatish Balay a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1663f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1664603f58a4SSatish Balay a->cowners = a->rowners + a->size + 2; 1665603f58a4SSatish Balay PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int)); 1666416022c9SBarry Smith ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 16672ee70a88SLois Curfman McInnes if (oldmat->colmap) { 16680452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1669416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1670416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1671416022c9SBarry Smith } else a->colmap = 0; 16723f41c07dSBarry Smith if (oldmat->garray) { 16733f41c07dSBarry Smith len = ((Mat_SeqAIJ *) (oldmat->B->data))->n; 16743f41c07dSBarry Smith a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray); 1675464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 16763f41c07dSBarry Smith if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1677416022c9SBarry Smith } else a->garray = 0; 1678d6dfbf8fSBarry Smith 1679416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1680416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1681a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1682416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 1683416022c9SBarry Smith ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 1684416022c9SBarry Smith PLogObjectParent(mat,a->A); 1685416022c9SBarry Smith ierr = MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr); 1686416022c9SBarry Smith PLogObjectParent(mat,a->B); 16875dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 168825cdf11fSBarry Smith if (flg) { 1689682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1690682d7d0cSBarry Smith } 16918a729477SBarry Smith *newmat = mat; 16928a729477SBarry Smith return 0; 16938a729477SBarry Smith } 1694416022c9SBarry Smith 169577c4ece6SBarry Smith #include "sys.h" 1696416022c9SBarry Smith 16975615d1e5SSatish Balay #undef __FUNC__ 16985615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ" 169919bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat) 1700416022c9SBarry Smith { 1701d65a2f8fSBarry Smith Mat A; 1702416022c9SBarry Smith int i, nz, ierr, j,rstart, rend, fd; 1703d65a2f8fSBarry Smith Scalar *vals,*svals; 170419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 1705416022c9SBarry Smith MPI_Status status; 170617699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1707d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 170819bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 1709416022c9SBarry Smith 171017699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 171117699dbbSLois Curfman McInnes if (!rank) { 171290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 171377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr); 1714e3372554SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object"); 1715416022c9SBarry Smith } 1716416022c9SBarry Smith 1717416022c9SBarry Smith MPI_Bcast(header+1,3,MPI_INT,0,comm); 1718416022c9SBarry Smith M = header[1]; N = header[2]; 1719416022c9SBarry Smith /* determine ownership of all rows */ 172017699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 17210452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1722416022c9SBarry Smith MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 1723416022c9SBarry Smith rowners[0] = 0; 172417699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1725416022c9SBarry Smith rowners[i] += rowners[i-1]; 1726416022c9SBarry Smith } 172717699dbbSLois Curfman McInnes rstart = rowners[rank]; 172817699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1729416022c9SBarry Smith 1730416022c9SBarry Smith /* distribute row lengths to all processors */ 17310452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1732416022c9SBarry Smith offlens = ourlens + (rend-rstart); 173317699dbbSLois Curfman McInnes if (!rank) { 17340452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 173577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 17360452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 173717699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1738d65a2f8fSBarry Smith MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 17390452661fSBarry Smith PetscFree(sndcounts); 1740416022c9SBarry Smith } 1741416022c9SBarry Smith else { 1742416022c9SBarry Smith MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 1743416022c9SBarry Smith } 1744416022c9SBarry Smith 174517699dbbSLois Curfman McInnes if (!rank) { 1746416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 17470452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1748cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 174917699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1750416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1751416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1752416022c9SBarry Smith } 1753416022c9SBarry Smith } 17540452661fSBarry Smith PetscFree(rowlengths); 1755416022c9SBarry Smith 1756416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1757416022c9SBarry Smith maxnz = 0; 175817699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 17590452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1760416022c9SBarry Smith } 17610452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1762416022c9SBarry Smith 1763416022c9SBarry Smith /* read in my part of the matrix column indices */ 1764416022c9SBarry Smith nz = procsnz[0]; 17650452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 176677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr); 1767d65a2f8fSBarry Smith 1768d65a2f8fSBarry Smith /* read in every one elses and ship off */ 176917699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1770d65a2f8fSBarry Smith nz = procsnz[i]; 177177c4ece6SBarry Smith ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 1772d65a2f8fSBarry Smith MPI_Send(cols,nz,MPI_INT,i,tag,comm); 1773d65a2f8fSBarry Smith } 17740452661fSBarry Smith PetscFree(cols); 1775416022c9SBarry Smith } 1776416022c9SBarry Smith else { 1777416022c9SBarry Smith /* determine buffer space needed for message */ 1778416022c9SBarry Smith nz = 0; 1779416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1780416022c9SBarry Smith nz += ourlens[i]; 1781416022c9SBarry Smith } 17820452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1783416022c9SBarry Smith 1784416022c9SBarry Smith /* receive message of column indices*/ 1785d65a2f8fSBarry Smith MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 1786416022c9SBarry Smith MPI_Get_count(&status,MPI_INT,&maxnz); 1787e3372554SBarry Smith if (maxnz != nz) SETERRQ(1,0,"something is wrong with file"); 1788416022c9SBarry Smith } 1789416022c9SBarry Smith 1790416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 1791cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 1792416022c9SBarry Smith jj = 0; 1793416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1794416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 1795d65a2f8fSBarry Smith if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1796416022c9SBarry Smith jj++; 1797416022c9SBarry Smith } 1798416022c9SBarry Smith } 1799d65a2f8fSBarry Smith 1800d65a2f8fSBarry Smith /* create our matrix */ 1801416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1802416022c9SBarry Smith ourlens[i] -= offlens[i]; 1803416022c9SBarry Smith } 1804d65a2f8fSBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 1805d65a2f8fSBarry Smith A = *newmat; 18066d4a8577SBarry Smith MatSetOption(A,MAT_COLUMNS_SORTED); 1807d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1808d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 1809d65a2f8fSBarry Smith } 1810416022c9SBarry Smith 181117699dbbSLois Curfman McInnes if (!rank) { 18120452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1813416022c9SBarry Smith 1814416022c9SBarry Smith /* read in my part of the matrix numerical values */ 1815416022c9SBarry Smith nz = procsnz[0]; 181677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 1817d65a2f8fSBarry Smith 1818d65a2f8fSBarry Smith /* insert into matrix */ 1819d65a2f8fSBarry Smith jj = rstart; 1820d65a2f8fSBarry Smith smycols = mycols; 1821d65a2f8fSBarry Smith svals = vals; 1822d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1823d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1824d65a2f8fSBarry Smith smycols += ourlens[i]; 1825d65a2f8fSBarry Smith svals += ourlens[i]; 1826d65a2f8fSBarry Smith jj++; 1827416022c9SBarry Smith } 1828416022c9SBarry Smith 1829d65a2f8fSBarry Smith /* read in other processors and ship out */ 183017699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1831416022c9SBarry Smith nz = procsnz[i]; 183277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 1833d65a2f8fSBarry Smith MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 1834416022c9SBarry Smith } 18350452661fSBarry Smith PetscFree(procsnz); 1836416022c9SBarry Smith } 1837d65a2f8fSBarry Smith else { 1838d65a2f8fSBarry Smith /* receive numeric values */ 18390452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1840416022c9SBarry Smith 1841d65a2f8fSBarry Smith /* receive message of values*/ 1842d65a2f8fSBarry Smith MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 1843d65a2f8fSBarry Smith MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 1844e3372554SBarry Smith if (maxnz != nz) SETERRQ(1,0,"something is wrong with file"); 1845d65a2f8fSBarry Smith 1846d65a2f8fSBarry Smith /* insert into matrix */ 1847d65a2f8fSBarry Smith jj = rstart; 1848d65a2f8fSBarry Smith smycols = mycols; 1849d65a2f8fSBarry Smith svals = vals; 1850d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1851d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1852d65a2f8fSBarry Smith smycols += ourlens[i]; 1853d65a2f8fSBarry Smith svals += ourlens[i]; 1854d65a2f8fSBarry Smith jj++; 1855d65a2f8fSBarry Smith } 1856d65a2f8fSBarry Smith } 18570452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1858d65a2f8fSBarry Smith 18596d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 18606d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1861416022c9SBarry Smith return 0; 1862416022c9SBarry Smith } 1863