1cb512458SBarry Smith #ifndef lint 2*f09e8eb9SSatish Balay static char vcid[] = "$Id: mpiaij.c,v 1.206 1997/05/09 22:34:48 balay 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 325615d1e5SSatish Balay #undef __FUNC__ 335615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_MPIAIJ" 348f6be9afSLois Curfman McInnes int MatGetRowIJ_MPIAIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int **ja, 353b2fbd54SBarry Smith PetscTruth *done) 36299609e3SLois Curfman McInnes { 37299609e3SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 38299609e3SLois Curfman McInnes int ierr; 3917699dbbSLois Curfman McInnes if (aij->size == 1) { 403b2fbd54SBarry Smith ierr = MatGetRowIJ(aij->A,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr); 41e3372554SBarry Smith } else SETERRQ(1,0,"not supported in parallel"); 423b2fbd54SBarry Smith return 0; 433b2fbd54SBarry Smith } 443b2fbd54SBarry Smith 455615d1e5SSatish Balay #undef __FUNC__ 465eea60f9SBarry Smith #define __FUNC__ "MatRestoreRowIJ_MPIAIJ" /* ADIC Ignore */ 478f6be9afSLois Curfman McInnes int MatRestoreRowIJ_MPIAIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int **ja, 483b2fbd54SBarry Smith PetscTruth *done) 493b2fbd54SBarry Smith { 503b2fbd54SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 513b2fbd54SBarry Smith int ierr; 523b2fbd54SBarry Smith if (aij->size == 1) { 533b2fbd54SBarry Smith ierr = MatRestoreRowIJ(aij->A,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr); 54e3372554SBarry Smith } else SETERRQ(1,0,"not supported in parallel"); 55299609e3SLois Curfman McInnes return 0; 56299609e3SLois Curfman McInnes } 57299609e3SLois Curfman McInnes 580520107fSSatish Balay #define CHUNKSIZE 15 5930770e4dSSatish Balay #define MatSetValues_SeqAIJ_A_Private(row,col,value,addv) \ 600520107fSSatish Balay { \ 610520107fSSatish Balay \ 620520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 6330770e4dSSatish Balay rmax = aimax[row]; nrow = ailen[row]; \ 64f5e9677aSSatish Balay col1 = col - shift; \ 65f5e9677aSSatish Balay \ 660520107fSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 67f5e9677aSSatish Balay if (rp[_i] > col1) break; \ 68f5e9677aSSatish Balay if (rp[_i] == col1) { \ 690520107fSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 700520107fSSatish Balay else ap[_i] = value; \ 7130770e4dSSatish Balay goto a_noinsert; \ 720520107fSSatish Balay } \ 730520107fSSatish Balay } \ 7489280ab3SLois Curfman McInnes if (nonew == 1) goto a_noinsert; \ 7511ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 760520107fSSatish Balay if (nrow >= rmax) { \ 770520107fSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 780520107fSSatish Balay int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \ 790520107fSSatish Balay Scalar *new_a; \ 800520107fSSatish Balay \ 8111ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 8289280ab3SLois Curfman McInnes \ 830520107fSSatish Balay /* malloc new storage space */ \ 840520107fSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \ 850520107fSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 860520107fSSatish Balay new_j = (int *) (new_a + new_nz); \ 870520107fSSatish Balay new_i = new_j + new_nz; \ 880520107fSSatish Balay \ 890520107fSSatish Balay /* copy over old data into new slots */ \ 900520107fSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \ 910520107fSSatish Balay for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \ 920520107fSSatish Balay PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \ 930520107fSSatish Balay len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \ 940520107fSSatish Balay PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \ 950520107fSSatish Balay len*sizeof(int)); \ 960520107fSSatish Balay PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \ 970520107fSSatish Balay PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \ 980520107fSSatish Balay len*sizeof(Scalar)); \ 990520107fSSatish Balay /* free up old matrix storage */ \ 100f5e9677aSSatish Balay \ 1010520107fSSatish Balay PetscFree(a->a); \ 1020520107fSSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \ 1030520107fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; \ 1040520107fSSatish Balay a->singlemalloc = 1; \ 1050520107fSSatish Balay \ 1060520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 10730770e4dSSatish Balay rmax = aimax[row] = aimax[row] + CHUNKSIZE; \ 1080520107fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 1090520107fSSatish Balay a->maxnz += CHUNKSIZE; \ 1100520107fSSatish Balay a->reallocs++; \ 1110520107fSSatish Balay } \ 1120520107fSSatish Balay N = nrow++ - 1; a->nz++; \ 1130520107fSSatish Balay /* shift up all the later entries in this row */ \ 1140520107fSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 1150520107fSSatish Balay rp[ii+1] = rp[ii]; \ 1160520107fSSatish Balay ap[ii+1] = ap[ii]; \ 1170520107fSSatish Balay } \ 118f5e9677aSSatish Balay rp[_i] = col1; \ 1190520107fSSatish Balay ap[_i] = value; \ 12030770e4dSSatish Balay a_noinsert: ; \ 1210520107fSSatish Balay ailen[row] = nrow; \ 1220520107fSSatish Balay } 1230a198c4cSBarry Smith 12430770e4dSSatish Balay #define MatSetValues_SeqAIJ_B_Private(row,col,value,addv) \ 12530770e4dSSatish Balay { \ 12630770e4dSSatish Balay \ 12730770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 12830770e4dSSatish Balay rmax = bimax[row]; nrow = bilen[row]; \ 12930770e4dSSatish Balay col1 = col - shift; \ 13030770e4dSSatish Balay \ 13130770e4dSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 13230770e4dSSatish Balay if (rp[_i] > col1) break; \ 13330770e4dSSatish Balay if (rp[_i] == col1) { \ 13430770e4dSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 13530770e4dSSatish Balay else ap[_i] = value; \ 13630770e4dSSatish Balay goto b_noinsert; \ 13730770e4dSSatish Balay } \ 13830770e4dSSatish Balay } \ 13989280ab3SLois Curfman McInnes if (nonew == 1) goto b_noinsert; \ 14011ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 14130770e4dSSatish Balay if (nrow >= rmax) { \ 14230770e4dSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 14374c639caSSatish Balay int new_nz = bi[b->m] + CHUNKSIZE,len,*new_i,*new_j; \ 14430770e4dSSatish Balay Scalar *new_a; \ 14530770e4dSSatish Balay \ 14611ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \ 14789280ab3SLois Curfman McInnes \ 14830770e4dSSatish Balay /* malloc new storage space */ \ 14974c639caSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(b->m+1)*sizeof(int); \ 15030770e4dSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 15130770e4dSSatish Balay new_j = (int *) (new_a + new_nz); \ 15230770e4dSSatish Balay new_i = new_j + new_nz; \ 15330770e4dSSatish Balay \ 15430770e4dSSatish Balay /* copy over old data into new slots */ \ 15530770e4dSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = bi[ii];} \ 15674c639caSSatish Balay for ( ii=row+1; ii<b->m+1; ii++ ) {new_i[ii] = bi[ii]+CHUNKSIZE;} \ 15730770e4dSSatish Balay PetscMemcpy(new_j,bj,(bi[row]+nrow+shift)*sizeof(int)); \ 15830770e4dSSatish Balay len = (new_nz - CHUNKSIZE - bi[row] - nrow - shift); \ 15930770e4dSSatish Balay PetscMemcpy(new_j+bi[row]+shift+nrow+CHUNKSIZE,bj+bi[row]+shift+nrow, \ 16030770e4dSSatish Balay len*sizeof(int)); \ 16130770e4dSSatish Balay PetscMemcpy(new_a,ba,(bi[row]+nrow+shift)*sizeof(Scalar)); \ 16230770e4dSSatish Balay PetscMemcpy(new_a+bi[row]+shift+nrow+CHUNKSIZE,ba+bi[row]+shift+nrow, \ 16330770e4dSSatish Balay len*sizeof(Scalar)); \ 16430770e4dSSatish Balay /* free up old matrix storage */ \ 16530770e4dSSatish Balay \ 16674c639caSSatish Balay PetscFree(b->a); \ 16774c639caSSatish Balay if (!b->singlemalloc) {PetscFree(b->i);PetscFree(b->j);} \ 16874c639caSSatish Balay ba = b->a = new_a; bi = b->i = new_i; bj = b->j = new_j; \ 16974c639caSSatish Balay b->singlemalloc = 1; \ 17030770e4dSSatish Balay \ 17130770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 17230770e4dSSatish Balay rmax = bimax[row] = bimax[row] + CHUNKSIZE; \ 17374c639caSSatish Balay PLogObjectMemory(B,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 17474c639caSSatish Balay b->maxnz += CHUNKSIZE; \ 17574c639caSSatish Balay b->reallocs++; \ 17630770e4dSSatish Balay } \ 17774c639caSSatish Balay N = nrow++ - 1; b->nz++; \ 17830770e4dSSatish Balay /* shift up all the later entries in this row */ \ 17930770e4dSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 18030770e4dSSatish Balay rp[ii+1] = rp[ii]; \ 18130770e4dSSatish Balay ap[ii+1] = ap[ii]; \ 18230770e4dSSatish Balay } \ 18330770e4dSSatish Balay rp[_i] = col1; \ 18430770e4dSSatish Balay ap[_i] = value; \ 18530770e4dSSatish Balay b_noinsert: ; \ 18630770e4dSSatish Balay bilen[row] = nrow; \ 18730770e4dSSatish Balay } 18830770e4dSSatish Balay 1890520107fSSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode); 1905615d1e5SSatish Balay #undef __FUNC__ 1915615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIAIJ" 1928f6be9afSLois Curfman McInnes int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 1938a729477SBarry Smith { 19444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1954b0e389bSBarry Smith Scalar value; 1961eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 1971eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 198905e6a2fSBarry Smith int roworiented = aij->roworiented; 1998a729477SBarry Smith 2000520107fSSatish Balay /* Some Variables required in the macro */ 2014ee7247eSSatish Balay Mat A = aij->A; 2024ee7247eSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 20330770e4dSSatish Balay int *aimax = a->imax, *ai = a->i, *ailen = a->ilen,*aj = a->j; 20430770e4dSSatish Balay Scalar *aa = a->a; 20530770e4dSSatish Balay 20630770e4dSSatish Balay Mat B = aij->B; 20730770e4dSSatish Balay Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data; 20830770e4dSSatish Balay int *bimax = b->imax, *bi = b->i, *bilen = b->ilen,*bj = b->j; 20930770e4dSSatish Balay Scalar *ba = b->a; 21030770e4dSSatish Balay 2114ee7247eSSatish Balay int *rp,ii,nrow,_i,rmax, N, col1; 21230770e4dSSatish Balay int nonew = a->nonew,shift = a->indexshift; 21330770e4dSSatish Balay Scalar *ap; 2144ee7247eSSatish Balay 2158a729477SBarry Smith for ( i=0; i<m; i++ ) { 2160a198c4cSBarry Smith #if defined(PETSC_BOPT_g) 217e3372554SBarry Smith if (im[i] < 0) SETERRQ(1,0,"Negative row"); 218e3372554SBarry Smith if (im[i] >= aij->M) SETERRQ(1,0,"Row too large"); 2190a198c4cSBarry Smith #endif 2204b0e389bSBarry Smith if (im[i] >= rstart && im[i] < rend) { 2214b0e389bSBarry Smith row = im[i] - rstart; 2221eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 2234b0e389bSBarry Smith if (in[j] >= cstart && in[j] < cend){ 2244b0e389bSBarry Smith col = in[j] - cstart; 2254b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 22630770e4dSSatish Balay MatSetValues_SeqAIJ_A_Private(row,col,value,addv); 2270520107fSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2281eb62cbbSBarry Smith } 2290a198c4cSBarry Smith #if defined(PETSC_BOPT_g) 230e3372554SBarry Smith else if (in[j] < 0) {SETERRQ(1,0,"Negative column");} 231e3372554SBarry Smith else if (in[j] >= aij->N) {SETERRQ(1,0,"Col too large");} 2320a198c4cSBarry Smith #endif 2331eb62cbbSBarry Smith else { 234227d817aSBarry Smith if (mat->was_assembled) { 235905e6a2fSBarry Smith if (!aij->colmap) { 236905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 237905e6a2fSBarry Smith } 238905e6a2fSBarry Smith col = aij->colmap[in[j]] - 1; 239ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 2402493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 2414b0e389bSBarry Smith col = in[j]; 2429bf004c3SSatish Balay /* Reinitialize the variables required by MatSetValues_SeqAIJ_B_Private() */ 243f9508a3cSSatish Balay B = aij->B; 244f9508a3cSSatish Balay b = (Mat_SeqAIJ *) B->data; 245f9508a3cSSatish Balay bimax = b->imax; bi = b->i; bilen = b->ilen; bj = b->j; 246f9508a3cSSatish Balay ba = b->a; 247d6dfbf8fSBarry Smith } 2489e25ed09SBarry Smith } 2494b0e389bSBarry Smith else col = in[j]; 2504b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 25130770e4dSSatish Balay MatSetValues_SeqAIJ_B_Private(row,col,value,addv); 25230770e4dSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2531eb62cbbSBarry Smith } 2541eb62cbbSBarry Smith } 2551eb62cbbSBarry Smith } 2561eb62cbbSBarry Smith else { 25790f02eecSBarry Smith if (roworiented && !aij->donotstash) { 2584b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr); 2594b0e389bSBarry Smith } 2604b0e389bSBarry Smith else { 26190f02eecSBarry Smith if (!aij->donotstash) { 2624b0e389bSBarry Smith row = im[i]; 2634b0e389bSBarry Smith for ( j=0; j<n; j++ ) { 2644b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr); 2654b0e389bSBarry Smith } 2664b0e389bSBarry Smith } 2671eb62cbbSBarry Smith } 2688a729477SBarry Smith } 26990f02eecSBarry Smith } 2708a729477SBarry Smith return 0; 2718a729477SBarry Smith } 2728a729477SBarry Smith 2735615d1e5SSatish Balay #undef __FUNC__ 2745615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIAIJ" 2758f6be9afSLois Curfman McInnes int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 276b49de8d1SLois Curfman McInnes { 277b49de8d1SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 278b49de8d1SLois Curfman McInnes int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 279b49de8d1SLois Curfman McInnes int cstart = aij->cstart, cend = aij->cend,row,col; 280b49de8d1SLois Curfman McInnes 281b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 282e3372554SBarry Smith if (idxm[i] < 0) SETERRQ(1,0,"Negative row"); 283e3372554SBarry Smith if (idxm[i] >= aij->M) SETERRQ(1,0,"Row too large"); 284b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 285b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 286b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 287e3372554SBarry Smith if (idxn[j] < 0) SETERRQ(1,0,"Negative column"); 288e3372554SBarry Smith if (idxn[j] >= aij->N) SETERRQ(1,0,"Col too large"); 289b49de8d1SLois Curfman McInnes if (idxn[j] >= cstart && idxn[j] < cend){ 290b49de8d1SLois Curfman McInnes col = idxn[j] - cstart; 291b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 292b49de8d1SLois Curfman McInnes } 293b49de8d1SLois Curfman McInnes else { 294905e6a2fSBarry Smith if (!aij->colmap) { 295905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 296905e6a2fSBarry Smith } 297905e6a2fSBarry Smith col = aij->colmap[idxn[j]] - 1; 298e60e1c95SSatish Balay if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0; 299d9d09a02SSatish Balay else { 300b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 301b49de8d1SLois Curfman McInnes } 302b49de8d1SLois Curfman McInnes } 303b49de8d1SLois Curfman McInnes } 304d9d09a02SSatish Balay } 305b49de8d1SLois Curfman McInnes else { 306e3372554SBarry Smith SETERRQ(1,0,"Only local values currently supported"); 307b49de8d1SLois Curfman McInnes } 308b49de8d1SLois Curfman McInnes } 309b49de8d1SLois Curfman McInnes return 0; 310b49de8d1SLois Curfman McInnes } 311b49de8d1SLois Curfman McInnes 3125615d1e5SSatish Balay #undef __FUNC__ 3135615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ" 3148f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode) 3158a729477SBarry Smith { 31644a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 317d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 31817699dbbSLois Curfman McInnes int size = aij->size, *owners = aij->rowners; 31917699dbbSLois Curfman McInnes int rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr; 3201eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3216abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 3221eb62cbbSBarry Smith InsertMode addv; 3231eb62cbbSBarry Smith Scalar *rvalues,*svalues; 3241eb62cbbSBarry Smith 3251eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 32647794344SBarry Smith MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm); 327dbb450caSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 328e3372554SBarry Smith SETERRQ(1,0,"Some processors inserted others added"); 3291eb62cbbSBarry Smith } 33047794344SBarry Smith mat->insertmode = addv; /* in case this processor had no cache */ 3311eb62cbbSBarry Smith 3321eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3330452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 334cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 3350452661fSBarry Smith owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 3361eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3371eb62cbbSBarry Smith idx = aij->stash.idx[i]; 33817699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 3391eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3401eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 3418a729477SBarry Smith } 3428a729477SBarry Smith } 3438a729477SBarry Smith } 34417699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 3451eb62cbbSBarry Smith 3461eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3470452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 34817699dbbSLois Curfman McInnes MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm); 34917699dbbSLois Curfman McInnes nreceives = work[rank]; 35017699dbbSLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 35117699dbbSLois Curfman McInnes nmax = work[rank]; 3520452661fSBarry Smith PetscFree(work); 3531eb62cbbSBarry Smith 3541eb62cbbSBarry Smith /* post receives: 3551eb62cbbSBarry Smith 1) each message will consist of ordered pairs 3561eb62cbbSBarry Smith (global index,value) we store the global index as a double 357d6dfbf8fSBarry Smith to simplify the message passing. 3581eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 3591eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 3601eb62cbbSBarry Smith this is a lot of wasted space. 3611eb62cbbSBarry Smith 3621eb62cbbSBarry Smith 3631eb62cbbSBarry Smith This could be done better. 3641eb62cbbSBarry Smith */ 3650452661fSBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 36678b31e54SBarry Smith CHKPTRQ(rvalues); 3670452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request)); 36878b31e54SBarry Smith CHKPTRQ(recv_waits); 3691eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 370416022c9SBarry Smith MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 3711eb62cbbSBarry Smith comm,recv_waits+i); 3721eb62cbbSBarry Smith } 3731eb62cbbSBarry Smith 3741eb62cbbSBarry Smith /* do sends: 3751eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3761eb62cbbSBarry Smith the ith processor 3771eb62cbbSBarry Smith */ 3780452661fSBarry Smith svalues = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 3790452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 38078b31e54SBarry Smith CHKPTRQ(send_waits); 3810452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 3821eb62cbbSBarry Smith starts[0] = 0; 38317699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3841eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 3851eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 3861eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 3871eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 3881eb62cbbSBarry Smith } 3890452661fSBarry Smith PetscFree(owner); 3901eb62cbbSBarry Smith starts[0] = 0; 39117699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3921eb62cbbSBarry Smith count = 0; 39317699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 3941eb62cbbSBarry Smith if (procs[i]) { 395416022c9SBarry Smith MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag, 3961eb62cbbSBarry Smith comm,send_waits+count++); 3971eb62cbbSBarry Smith } 3981eb62cbbSBarry Smith } 3990452661fSBarry Smith PetscFree(starts); PetscFree(nprocs); 4001eb62cbbSBarry Smith 4011eb62cbbSBarry Smith /* Free cache space */ 40255908cccSLois Curfman McInnes PLogInfo(mat,"MatAssemblyBegin_MPIAIJ:Number of off-processor values %d\n",aij->stash.n); 40378b31e54SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 4041eb62cbbSBarry Smith 4051eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 4061eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 4071eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 4081eb62cbbSBarry Smith aij->rmax = nmax; 4091eb62cbbSBarry Smith 4101eb62cbbSBarry Smith return 0; 4111eb62cbbSBarry Smith } 41244a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat); 4131eb62cbbSBarry Smith 4145615d1e5SSatish Balay #undef __FUNC__ 4155615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ" 4168f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode) 4171eb62cbbSBarry Smith { 41844a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 4191eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 420416022c9SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr; 421905e6a2fSBarry Smith int row,col,other_disassembled; 4221eb62cbbSBarry Smith Scalar *values,val; 42347794344SBarry Smith InsertMode addv = mat->insertmode; 4241eb62cbbSBarry Smith 4251eb62cbbSBarry Smith /* wait on receives */ 4261eb62cbbSBarry Smith while (count) { 427d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status); 4281eb62cbbSBarry Smith /* unpack receives into our local space */ 429d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 430c60448a5SBarry Smith MPI_Get_count(&recv_status,MPIU_SCALAR,&n); 4311eb62cbbSBarry Smith n = n/3; 4321eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 433227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - aij->rstart; 434227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 4351eb62cbbSBarry Smith val = values[3*i+2]; 4361eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 4371eb62cbbSBarry Smith col -= aij->cstart; 4381eb62cbbSBarry Smith MatSetValues(aij->A,1,&row,1,&col,&val,addv); 4391eb62cbbSBarry Smith } 4401eb62cbbSBarry Smith else { 44155a1b374SBarry Smith if (mat->was_assembled || mat->assembled) { 442905e6a2fSBarry Smith if (!aij->colmap) { 443905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 444905e6a2fSBarry Smith } 445905e6a2fSBarry Smith col = aij->colmap[col] - 1; 446ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 4472493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 448227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 449d6dfbf8fSBarry Smith } 4509e25ed09SBarry Smith } 4511eb62cbbSBarry Smith MatSetValues(aij->B,1,&row,1,&col,&val,addv); 4521eb62cbbSBarry Smith } 4531eb62cbbSBarry Smith } 4541eb62cbbSBarry Smith count--; 4551eb62cbbSBarry Smith } 4560452661fSBarry Smith PetscFree(aij->recv_waits); PetscFree(aij->rvalues); 4571eb62cbbSBarry Smith 4581eb62cbbSBarry Smith /* wait on sends */ 4591eb62cbbSBarry Smith if (aij->nsends) { 4600a198c4cSBarry Smith send_status = (MPI_Status *) PetscMalloc(aij->nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 4611eb62cbbSBarry Smith MPI_Waitall(aij->nsends,aij->send_waits,send_status); 4620452661fSBarry Smith PetscFree(send_status); 4631eb62cbbSBarry Smith } 4640452661fSBarry Smith PetscFree(aij->send_waits); PetscFree(aij->svalues); 4651eb62cbbSBarry Smith 46678b31e54SBarry Smith ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr); 46778b31e54SBarry Smith ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr); 4681eb62cbbSBarry Smith 4692493cbb0SBarry Smith /* determine if any processor has disassembled, if so we must 4702493cbb0SBarry Smith also disassemble ourselfs, in order that we may reassemble. */ 471227d817aSBarry Smith MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm); 472227d817aSBarry Smith if (mat->was_assembled && !other_disassembled) { 4732493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 4742493cbb0SBarry Smith } 4752493cbb0SBarry Smith 4766d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 47778b31e54SBarry Smith ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr); 4785e42470aSBarry Smith } 47978b31e54SBarry Smith ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr); 48078b31e54SBarry Smith ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr); 4815e42470aSBarry Smith 4827a0afa10SBarry Smith if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;} 4838a729477SBarry Smith return 0; 4848a729477SBarry Smith } 4858a729477SBarry Smith 4865615d1e5SSatish Balay #undef __FUNC__ 4875615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ" 4888f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A) 4891eb62cbbSBarry Smith { 49044a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 491dbd7a890SLois Curfman McInnes int ierr; 49278b31e54SBarry Smith ierr = MatZeroEntries(l->A); CHKERRQ(ierr); 49378b31e54SBarry Smith ierr = MatZeroEntries(l->B); CHKERRQ(ierr); 4941eb62cbbSBarry Smith return 0; 4951eb62cbbSBarry Smith } 4961eb62cbbSBarry Smith 4971eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 4981eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 4991eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 5001eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 5011eb62cbbSBarry Smith routine. 5021eb62cbbSBarry Smith */ 5035615d1e5SSatish Balay #undef __FUNC__ 5045615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ" 5058f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 5061eb62cbbSBarry Smith { 50744a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 50817699dbbSLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 5096abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 51017699dbbSLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 5115392566eSBarry Smith int *rvalues,tag = A->tag,count,base,slen,n,*source; 512d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 513d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 5141eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 5151eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 5161eb62cbbSBarry Smith IS istmp; 5171eb62cbbSBarry Smith 51877c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 51978b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5201eb62cbbSBarry Smith 5211eb62cbbSBarry Smith /* first count number of contributors to each processor */ 5220452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 523cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 5240452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 5251eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5261eb62cbbSBarry Smith idx = rows[i]; 5271eb62cbbSBarry Smith found = 0; 52817699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 5291eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 5301eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 5311eb62cbbSBarry Smith } 5321eb62cbbSBarry Smith } 533e3372554SBarry Smith if (!found) SETERRQ(1,0,"Index out of range"); 5341eb62cbbSBarry Smith } 53517699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 5361eb62cbbSBarry Smith 5371eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 5380452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 53917699dbbSLois Curfman McInnes MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm); 54017699dbbSLois Curfman McInnes nrecvs = work[rank]; 54117699dbbSLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 54217699dbbSLois Curfman McInnes nmax = work[rank]; 5430452661fSBarry Smith PetscFree(work); 5441eb62cbbSBarry Smith 5451eb62cbbSBarry Smith /* post receives: */ 5460452661fSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 54778b31e54SBarry Smith CHKPTRQ(rvalues); 5480452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request)); 54978b31e54SBarry Smith CHKPTRQ(recv_waits); 5501eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 551416022c9SBarry Smith MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i); 5521eb62cbbSBarry Smith } 5531eb62cbbSBarry Smith 5541eb62cbbSBarry Smith /* do sends: 5551eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 5561eb62cbbSBarry Smith the ith processor 5571eb62cbbSBarry Smith */ 5580452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 5590452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 56078b31e54SBarry Smith CHKPTRQ(send_waits); 5610452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 5621eb62cbbSBarry Smith starts[0] = 0; 56317699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5641eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 5651eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 5661eb62cbbSBarry Smith } 5671eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 5681eb62cbbSBarry Smith 5691eb62cbbSBarry Smith starts[0] = 0; 57017699dbbSLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 5711eb62cbbSBarry Smith count = 0; 57217699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 5731eb62cbbSBarry Smith if (procs[i]) { 574416022c9SBarry Smith MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++); 5751eb62cbbSBarry Smith } 5761eb62cbbSBarry Smith } 5770452661fSBarry Smith PetscFree(starts); 5781eb62cbbSBarry Smith 57917699dbbSLois Curfman McInnes base = owners[rank]; 5801eb62cbbSBarry Smith 5811eb62cbbSBarry Smith /* wait on receives */ 5820452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 5831eb62cbbSBarry Smith source = lens + nrecvs; 5841eb62cbbSBarry Smith count = nrecvs; slen = 0; 5851eb62cbbSBarry Smith while (count) { 586d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 5871eb62cbbSBarry Smith /* unpack receives into our local space */ 5881eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_INT,&n); 589d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 590d6dfbf8fSBarry Smith lens[imdex] = n; 5911eb62cbbSBarry Smith slen += n; 5921eb62cbbSBarry Smith count--; 5931eb62cbbSBarry Smith } 5940452661fSBarry Smith PetscFree(recv_waits); 5951eb62cbbSBarry Smith 5961eb62cbbSBarry Smith /* move the data into the send scatter */ 5970452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 5981eb62cbbSBarry Smith count = 0; 5991eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 6001eb62cbbSBarry Smith values = rvalues + i*nmax; 6011eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 6021eb62cbbSBarry Smith lrows[count++] = values[j] - base; 6031eb62cbbSBarry Smith } 6041eb62cbbSBarry Smith } 6050452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 6060452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 6071eb62cbbSBarry Smith 6081eb62cbbSBarry Smith /* actually zap the local rows */ 609029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 610464493b3SBarry Smith PLogObjectParent(A,istmp); 6110452661fSBarry Smith PetscFree(lrows); 61278b31e54SBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 61378b31e54SBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr); 61478b31e54SBarry Smith ierr = ISDestroy(istmp); CHKERRQ(ierr); 6151eb62cbbSBarry Smith 6161eb62cbbSBarry Smith /* wait on sends */ 6171eb62cbbSBarry Smith if (nsends) { 6180452661fSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status)); 61978b31e54SBarry Smith CHKPTRQ(send_status); 6201eb62cbbSBarry Smith MPI_Waitall(nsends,send_waits,send_status); 6210452661fSBarry Smith PetscFree(send_status); 6221eb62cbbSBarry Smith } 6230452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 6241eb62cbbSBarry Smith 6251eb62cbbSBarry Smith return 0; 6261eb62cbbSBarry Smith } 6271eb62cbbSBarry Smith 6285615d1e5SSatish Balay #undef __FUNC__ 6295615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ" 6308f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy) 6311eb62cbbSBarry Smith { 632416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 633fbd6ef76SBarry Smith int ierr,nt; 634416022c9SBarry Smith 635a2ce50c7SBarry Smith ierr = VecGetLocalSize(xx,&nt); CHKERRQ(ierr); 636fbd6ef76SBarry Smith if (nt != a->n) { 637f40265daSLois Curfman McInnes SETERRQ(1,0,"Incompatible partition of A and xx"); 638fbd6ef76SBarry Smith } 63943a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 64038597bd4SSatish Balay ierr = (*a->A->ops.mult)(a->A,xx,yy); CHKERRQ(ierr); 64143a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 64238597bd4SSatish Balay ierr = (*a->B->ops.multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr); 6431eb62cbbSBarry Smith return 0; 6441eb62cbbSBarry Smith } 6451eb62cbbSBarry Smith 6465615d1e5SSatish Balay #undef __FUNC__ 6475615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ" 6488f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 649da3a660dSBarry Smith { 650416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 651da3a660dSBarry Smith int ierr; 65243a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 65338597bd4SSatish Balay ierr = (*a->A->ops.multadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 65443a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 65538597bd4SSatish Balay ierr = (*a->B->ops.multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr); 656da3a660dSBarry Smith return 0; 657da3a660dSBarry Smith } 658da3a660dSBarry Smith 6595615d1e5SSatish Balay #undef __FUNC__ 6605615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ" 6618f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy) 662da3a660dSBarry Smith { 663416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 664da3a660dSBarry Smith int ierr; 665da3a660dSBarry Smith 666da3a660dSBarry Smith /* do nondiagonal part */ 66738597bd4SSatish Balay ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 668da3a660dSBarry Smith /* send it on its way */ 669537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 670da3a660dSBarry Smith /* do local part */ 67138597bd4SSatish Balay ierr = (*a->A->ops.multtrans)(a->A,xx,yy); CHKERRQ(ierr); 672da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 673da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 674da3a660dSBarry Smith /* but is not perhaps always true. */ 675537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 676da3a660dSBarry Smith return 0; 677da3a660dSBarry Smith } 678da3a660dSBarry Smith 6795615d1e5SSatish Balay #undef __FUNC__ 6805615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ" 6818f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 682da3a660dSBarry Smith { 683416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 684da3a660dSBarry Smith int ierr; 685da3a660dSBarry Smith 686da3a660dSBarry Smith /* do nondiagonal part */ 68738597bd4SSatish Balay ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 688da3a660dSBarry Smith /* send it on its way */ 689537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 690da3a660dSBarry Smith /* do local part */ 69138597bd4SSatish Balay ierr = (*a->A->ops.multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 692da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 693da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 694da3a660dSBarry Smith /* but is not perhaps always true. */ 6950a198c4cSBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 696da3a660dSBarry Smith return 0; 697da3a660dSBarry Smith } 698da3a660dSBarry Smith 6991eb62cbbSBarry Smith /* 7001eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 7011eb62cbbSBarry Smith diagonal block 7021eb62cbbSBarry Smith */ 7035615d1e5SSatish Balay #undef __FUNC__ 7045615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ" 7058f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v) 7061eb62cbbSBarry Smith { 707416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 70844cd7ae7SLois Curfman McInnes if (a->M != a->N) 709e3372554SBarry Smith SETERRQ(1,0,"Supports only square matrix where A->A is diag block"); 7105baf8537SBarry Smith if (a->rstart != a->cstart || a->rend != a->cend) { 711e3372554SBarry Smith SETERRQ(1,0,"row partition must equal col partition"); } 712416022c9SBarry Smith return MatGetDiagonal(a->A,v); 7131eb62cbbSBarry Smith } 7141eb62cbbSBarry Smith 7155615d1e5SSatish Balay #undef __FUNC__ 7165615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ" 7178f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A) 718052efed2SBarry Smith { 719052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 720052efed2SBarry Smith int ierr; 721052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 722052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 723052efed2SBarry Smith return 0; 724052efed2SBarry Smith } 725052efed2SBarry Smith 7265615d1e5SSatish Balay #undef __FUNC__ 7275eea60f9SBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" /* ADIC Ignore */ 7288f6be9afSLois Curfman McInnes int MatDestroy_MPIAIJ(PetscObject obj) 7291eb62cbbSBarry Smith { 7301eb62cbbSBarry Smith Mat mat = (Mat) obj; 73144a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 7321eb62cbbSBarry Smith int ierr; 733a5a9c739SBarry Smith #if defined(PETSC_LOG) 7346e35fa57SLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N); 735a5a9c739SBarry Smith #endif 7360452661fSBarry Smith PetscFree(aij->rowners); 73778b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 73878b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 7390452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 7400452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 7411eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 742a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 7437a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 7440452661fSBarry Smith PetscFree(aij); 74590f02eecSBarry Smith if (mat->mapping) { 74690f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 74790f02eecSBarry Smith } 748a5a9c739SBarry Smith PLogObjectDestroy(mat); 7490452661fSBarry Smith PetscHeaderDestroy(mat); 7501eb62cbbSBarry Smith return 0; 7511eb62cbbSBarry Smith } 752ee50ffe9SBarry Smith 7535615d1e5SSatish Balay #undef __FUNC__ 7545eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" /* ADIC Ignore */ 7558f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 7561eb62cbbSBarry Smith { 757416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 758416022c9SBarry Smith int ierr; 759416022c9SBarry Smith 76017699dbbSLois Curfman McInnes if (aij->size == 1) { 761416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 762416022c9SBarry Smith } 763e3372554SBarry Smith else SETERRQ(1,0,"Only uniprocessor output supported"); 764416022c9SBarry Smith return 0; 765416022c9SBarry Smith } 766416022c9SBarry Smith 7675615d1e5SSatish Balay #undef __FUNC__ 7685eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab" /* ADIC Ignore */ 7698f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer) 770416022c9SBarry Smith { 77144a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 772dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 773a56f8943SBarry Smith int ierr, format,shift = C->indexshift,rank; 774d636dbe3SBarry Smith FILE *fd; 77519bcc07fSBarry Smith ViewerType vtype; 776416022c9SBarry Smith 77719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 77819bcc07fSBarry Smith if (vtype == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) { 77990ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 7800a198c4cSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 7814e220ebcSLois Curfman McInnes MatInfo info; 7824e220ebcSLois Curfman McInnes int flg; 783a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 78490ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 7854e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 78695e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 78777c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 78895e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 7894e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 79095e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 7914e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 7924e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->A,MAT_LOCAL,&info); 7934e220ebcSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used); 7944e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->B,MAT_LOCAL,&info); 7954e220ebcSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used); 796a56f8943SBarry Smith fflush(fd); 79777c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 798a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 799416022c9SBarry Smith return 0; 800416022c9SBarry Smith } 8010a198c4cSBarry Smith else if (format == VIEWER_FORMAT_ASCII_INFO) { 80208480c60SBarry Smith return 0; 80308480c60SBarry Smith } 804416022c9SBarry Smith } 805416022c9SBarry Smith 80619bcc07fSBarry Smith if (vtype == DRAW_VIEWER) { 80719bcc07fSBarry Smith Draw draw; 80819bcc07fSBarry Smith PetscTruth isnull; 80919bcc07fSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 81019bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 81119bcc07fSBarry Smith } 81219bcc07fSBarry Smith 81319bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 81490ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 81577c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 816d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 81717699dbbSLois Curfman McInnes aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 8181eb62cbbSBarry Smith aij->cend); 81978b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 82078b31e54SBarry Smith ierr = MatView(aij->B,viewer); CHKERRQ(ierr); 821d13ab20cSBarry Smith fflush(fd); 82277c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 823d13ab20cSBarry Smith } 824416022c9SBarry Smith else { 825a56f8943SBarry Smith int size = aij->size; 826a56f8943SBarry Smith rank = aij->rank; 82717699dbbSLois Curfman McInnes if (size == 1) { 82878b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 82995373324SBarry Smith } 83095373324SBarry Smith else { 83195373324SBarry Smith /* assemble the entire matrix onto first processor. */ 83295373324SBarry Smith Mat A; 833ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 8342eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 83595373324SBarry Smith Scalar *a; 8362ee70a88SLois Curfman McInnes 83717699dbbSLois Curfman McInnes if (!rank) { 838b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 839c451ab25SLois Curfman McInnes CHKERRQ(ierr); 84095373324SBarry Smith } 84195373324SBarry Smith else { 842b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 843c451ab25SLois Curfman McInnes CHKERRQ(ierr); 84495373324SBarry Smith } 845464493b3SBarry Smith PLogObjectParent(mat,A); 846416022c9SBarry Smith 84795373324SBarry Smith /* copy over the A part */ 848ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 8492ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 85095373324SBarry Smith row = aij->rstart; 851dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 85295373324SBarry Smith for ( i=0; i<m; i++ ) { 853416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 85495373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 85595373324SBarry Smith } 8562ee70a88SLois Curfman McInnes aj = Aloc->j; 857dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 85895373324SBarry Smith 85995373324SBarry Smith /* copy over the B part */ 860ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 8612ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 86295373324SBarry Smith row = aij->rstart; 8630452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 864dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 86595373324SBarry Smith for ( i=0; i<m; i++ ) { 866416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 86795373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 86895373324SBarry Smith } 8690452661fSBarry Smith PetscFree(ct); 8706d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8716d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 87217699dbbSLois Curfman McInnes if (!rank) { 87378b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 87495373324SBarry Smith } 87578b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 87695373324SBarry Smith } 87795373324SBarry Smith } 8781eb62cbbSBarry Smith return 0; 8791eb62cbbSBarry Smith } 8801eb62cbbSBarry Smith 8815615d1e5SSatish Balay #undef __FUNC__ 8825eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ" /* ADIC Ignore */ 8838f6be9afSLois Curfman McInnes int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 884416022c9SBarry Smith { 885416022c9SBarry Smith Mat mat = (Mat) obj; 886416022c9SBarry Smith int ierr; 88719bcc07fSBarry Smith ViewerType vtype; 888416022c9SBarry Smith 88919bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 89019bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER || 89119bcc07fSBarry Smith vtype == DRAW_VIEWER || vtype == MATLAB_VIEWER) { 892d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 893416022c9SBarry Smith } 89419bcc07fSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 895416022c9SBarry Smith return MatView_MPIAIJ_Binary(mat,viewer); 896416022c9SBarry Smith } 897416022c9SBarry Smith return 0; 898416022c9SBarry Smith } 899416022c9SBarry Smith 9001eb62cbbSBarry Smith /* 9011eb62cbbSBarry Smith This has to provide several versions. 9021eb62cbbSBarry Smith 9031eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 9041eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 905d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 9061eb62cbbSBarry Smith */ 9075615d1e5SSatish Balay #undef __FUNC__ 9085615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ" 9098f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 910dbb450caSBarry Smith double fshift,int its,Vec xx) 9118a729477SBarry Smith { 91244a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 913d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 914ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 915c16cb8f2SBarry Smith Scalar *b,*x,*xs,*ls,d,*v,sum; 9166abc6512SBarry Smith int ierr,*idx, *diag; 917416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 9188a729477SBarry Smith 919d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 920dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 921dbb450caSBarry Smith ls = ls + shift; 922ec8511deSBarry Smith if (!A->diag) {if ((ierr = MatMarkDiag_SeqAIJ(AA))) return ierr;} 923d6dfbf8fSBarry Smith diag = A->diag; 924c16cb8f2SBarry Smith if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 925da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 92638597bd4SSatish Balay return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx); 927da3a660dSBarry Smith } 92843a90d84SBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 92978b31e54SBarry Smith CHKERRQ(ierr); 93043a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 93178b31e54SBarry Smith CHKERRQ(ierr); 932d6dfbf8fSBarry Smith while (its--) { 933d6dfbf8fSBarry Smith /* go down through the rows */ 934d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 935d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 936dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 937dbb450caSBarry Smith v = A->a + A->i[i] + shift; 938d6dfbf8fSBarry Smith sum = b[i]; 939d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 940dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 941d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 942dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 943dbb450caSBarry Smith v = B->a + B->i[i] + shift; 944d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 94555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 946d6dfbf8fSBarry Smith } 947d6dfbf8fSBarry Smith /* come up through the rows */ 948d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 949d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 950dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 951dbb450caSBarry Smith v = A->a + A->i[i] + shift; 952d6dfbf8fSBarry Smith sum = b[i]; 953d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 954dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 955d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 956dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 957dbb450caSBarry Smith v = B->a + B->i[i] + shift; 958d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 95955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 960d6dfbf8fSBarry Smith } 961d6dfbf8fSBarry Smith } 962d6dfbf8fSBarry Smith } 963d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 964da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 96538597bd4SSatish Balay return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx); 966da3a660dSBarry Smith } 96743a90d84SBarry Smith ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 96878b31e54SBarry Smith CHKERRQ(ierr); 96943a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx); 97078b31e54SBarry Smith CHKERRQ(ierr); 971d6dfbf8fSBarry Smith while (its--) { 972d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 973d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 974dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 975dbb450caSBarry Smith v = A->a + A->i[i] + shift; 976d6dfbf8fSBarry Smith sum = b[i]; 977d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 978dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 979d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 980dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 981dbb450caSBarry Smith v = B->a + B->i[i] + shift; 982d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 98355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 984d6dfbf8fSBarry Smith } 985d6dfbf8fSBarry Smith } 986d6dfbf8fSBarry Smith } 987d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 988da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 98938597bd4SSatish Balay return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx); 990da3a660dSBarry Smith } 99143a90d84SBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 99278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 99343a90d84SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD, 99478b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 995d6dfbf8fSBarry Smith while (its--) { 996d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 997d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 998dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 999dbb450caSBarry Smith v = A->a + A->i[i] + shift; 1000d6dfbf8fSBarry Smith sum = b[i]; 1001d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 1002dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 1003d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 1004dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 1005dbb450caSBarry Smith v = B->a + B->i[i] + shift; 1006d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 100755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 1008d6dfbf8fSBarry Smith } 1009d6dfbf8fSBarry Smith } 1010d6dfbf8fSBarry Smith } 1011c16cb8f2SBarry Smith else { 1012e3372554SBarry Smith SETERRQ(1,0,"Option not supported"); 1013c16cb8f2SBarry Smith } 10148a729477SBarry Smith return 0; 10158a729477SBarry Smith } 1016a66be287SLois Curfman McInnes 10175615d1e5SSatish Balay #undef __FUNC__ 10185eea60f9SBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" /* ADIC Ignore */ 10198f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info) 1020a66be287SLois Curfman McInnes { 1021a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1022a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 10234e220ebcSLois Curfman McInnes int ierr; 10244e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 1025a66be287SLois Curfman McInnes 10264e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 10274e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 10284e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 10294e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 10304e220ebcSLois Curfman McInnes info->block_size = 1.0; 10314e220ebcSLois Curfman McInnes ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr); 10324e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 10334e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 10344e220ebcSLois Curfman McInnes ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr); 10354e220ebcSLois Curfman McInnes isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded; 10364e220ebcSLois Curfman McInnes isend[3] += info->memory; isend[4] += info->mallocs; 1037a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 10384e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 10394e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 10404e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 10414e220ebcSLois Curfman McInnes info->memory = isend[3]; 10424e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 1043a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 10444e220ebcSLois Curfman McInnes MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm); 10454e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10464e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10474e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10484e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10494e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1050a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 10514e220ebcSLois Curfman McInnes MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm); 10524e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10534e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10544e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10554e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10564e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1057a66be287SLois Curfman McInnes } 10584e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 10594e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10604e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10614e220ebcSLois Curfman McInnes 1062a66be287SLois Curfman McInnes return 0; 1063a66be287SLois Curfman McInnes } 1064a66be287SLois Curfman McInnes 1065299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*); 1066299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*); 1067299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double); 1068299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *); 1069299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec); 1070299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec); 1071299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec); 1072299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec); 1073299609e3SLois Curfman McInnes 10745615d1e5SSatish Balay #undef __FUNC__ 10755eea60f9SBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" /* ADIC Ignore */ 10768f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op) 1077c74985f6SBarry Smith { 1078c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1079c74985f6SBarry Smith 10806d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 10816d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 10826da5968aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 1083c2653b3dSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 108496854ed6SLois Curfman McInnes op == MAT_NEW_NONZERO_ALLOCATION_ERROR || 1085c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR) { 1086b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 1087b1fbbac0SLois Curfman McInnes MatSetOption(a->B,op); 1088b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 1089aeafbbfcSLois Curfman McInnes a->roworiented = 1; 1090c0bbcb79SLois Curfman McInnes MatSetOption(a->A,op); 1091c0bbcb79SLois Curfman McInnes MatSetOption(a->B,op); 1092b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 10936da5968aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 10946d4a8577SBarry Smith op == MAT_SYMMETRIC || 10956d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 10966d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 109794a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_MPIAIJ:Option ignored\n"); 10986d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) { 10994b0e389bSBarry Smith a->roworiented = 0; 11004b0e389bSBarry Smith MatSetOption(a->A,op); 11014b0e389bSBarry Smith MatSetOption(a->B,op); 11022b362799SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 110390f02eecSBarry Smith a->donotstash = 1; 110490f02eecSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) 1105e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");} 1106c0bbcb79SLois Curfman McInnes else 1107e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"unknown option");} 1108c74985f6SBarry Smith return 0; 1109c74985f6SBarry Smith } 1110c74985f6SBarry Smith 11115615d1e5SSatish Balay #undef __FUNC__ 11125eea60f9SBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" /* ADIC Ignore */ 11138f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1114c74985f6SBarry Smith { 111544a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1116c74985f6SBarry Smith *m = mat->M; *n = mat->N; 1117c74985f6SBarry Smith return 0; 1118c74985f6SBarry Smith } 1119c74985f6SBarry Smith 11205615d1e5SSatish Balay #undef __FUNC__ 11215eea60f9SBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" /* ADIC Ignore */ 11228f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1123c74985f6SBarry Smith { 112444a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1125b7c46309SBarry Smith *m = mat->m; *n = mat->N; 1126c74985f6SBarry Smith return 0; 1127c74985f6SBarry Smith } 1128c74985f6SBarry Smith 11295615d1e5SSatish Balay #undef __FUNC__ 11305eea60f9SBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" /* ADIC Ignore */ 11318f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1132c74985f6SBarry Smith { 113344a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1134c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 1135c74985f6SBarry Smith return 0; 1136c74985f6SBarry Smith } 1137c74985f6SBarry Smith 11386d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11396d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11406d84be18SBarry Smith 11415615d1e5SSatish Balay #undef __FUNC__ 11425615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ" 11436d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 114439e00950SLois Curfman McInnes { 1145154123eaSLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 114670f0671dSBarry Smith Scalar *vworkA, *vworkB, **pvA, **pvB,*v_p; 1147154123eaSLois Curfman McInnes int i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart; 1148154123eaSLois Curfman McInnes int nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend; 114970f0671dSBarry Smith int *cmap, *idx_p; 115039e00950SLois Curfman McInnes 1151e3372554SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,0,"Already active"); 11527a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11537a0afa10SBarry Smith 115470f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11557a0afa10SBarry Smith /* 11567a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11577a0afa10SBarry Smith */ 11587a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 1159c16cb8f2SBarry Smith int max = 1,m = mat->m,tmp; 1160c16cb8f2SBarry Smith for ( i=0; i<m; i++ ) { 11617a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11627a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11637a0afa10SBarry Smith } 11647a0afa10SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar))); 11657a0afa10SBarry Smith CHKPTRQ(mat->rowvalues); 11667a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11677a0afa10SBarry Smith } 11687a0afa10SBarry Smith 1169e3372554SBarry Smith if (row < rstart || row >= rend) SETERRQ(1,0,"Only local rows") 1170abc0e9e4SLois Curfman McInnes lrow = row - rstart; 117139e00950SLois Curfman McInnes 1172154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1173154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1174154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 117538597bd4SSatish Balay ierr = (*mat->A->ops.getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 117638597bd4SSatish Balay ierr = (*mat->B->ops.getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1177154123eaSLois Curfman McInnes nztot = nzA + nzB; 1178154123eaSLois Curfman McInnes 117970f0671dSBarry Smith cmap = mat->garray; 1180154123eaSLois Curfman McInnes if (v || idx) { 1181154123eaSLois Curfman McInnes if (nztot) { 1182154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 118370f0671dSBarry Smith int imark = -1; 1184154123eaSLois Curfman McInnes if (v) { 118570f0671dSBarry Smith *v = v_p = mat->rowvalues; 118639e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 118770f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1188154123eaSLois Curfman McInnes else break; 1189154123eaSLois Curfman McInnes } 1190154123eaSLois Curfman McInnes imark = i; 119170f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 119270f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1193154123eaSLois Curfman McInnes } 1194154123eaSLois Curfman McInnes if (idx) { 119570f0671dSBarry Smith *idx = idx_p = mat->rowindices; 119670f0671dSBarry Smith if (imark > -1) { 119770f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 119870f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 119970f0671dSBarry Smith } 120070f0671dSBarry Smith } else { 1201154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 120270f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1203154123eaSLois Curfman McInnes else break; 1204154123eaSLois Curfman McInnes } 1205154123eaSLois Curfman McInnes imark = i; 120670f0671dSBarry Smith } 120770f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 120870f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 120939e00950SLois Curfman McInnes } 121039e00950SLois Curfman McInnes } 12111ca473b0SSatish Balay else { 12121ca473b0SSatish Balay if (idx) *idx = 0; 12131ca473b0SSatish Balay if (v) *v = 0; 12141ca473b0SSatish Balay } 1215154123eaSLois Curfman McInnes } 121639e00950SLois Curfman McInnes *nz = nztot; 121738597bd4SSatish Balay ierr = (*mat->A->ops.restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 121838597bd4SSatish Balay ierr = (*mat->B->ops.restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 121939e00950SLois Curfman McInnes return 0; 122039e00950SLois Curfman McInnes } 122139e00950SLois Curfman McInnes 12225615d1e5SSatish Balay #undef __FUNC__ 12235eea60f9SBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" /* ADIC Ignore */ 12246d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 122539e00950SLois Curfman McInnes { 12267a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 12277a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 1228e3372554SBarry Smith SETERRQ(1,0,"MatGetRow not called"); 12297a0afa10SBarry Smith } 12307a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 123139e00950SLois Curfman McInnes return 0; 123239e00950SLois Curfman McInnes } 123339e00950SLois Curfman McInnes 12345615d1e5SSatish Balay #undef __FUNC__ 12355615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ" 12368f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1237855ac2c5SLois Curfman McInnes { 1238855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1239ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1240416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1241416022c9SBarry Smith double sum = 0.0; 124204ca555eSLois Curfman McInnes Scalar *v; 124304ca555eSLois Curfman McInnes 124417699dbbSLois Curfman McInnes if (aij->size == 1) { 124514183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 124637fa93a5SLois Curfman McInnes } else { 124704ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 124804ca555eSLois Curfman McInnes v = amat->a; 124904ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 125004ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 125104ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 125204ca555eSLois Curfman McInnes #else 125304ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 125404ca555eSLois Curfman McInnes #endif 125504ca555eSLois Curfman McInnes } 125604ca555eSLois Curfman McInnes v = bmat->a; 125704ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 125804ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 125904ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 126004ca555eSLois Curfman McInnes #else 126104ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 126204ca555eSLois Curfman McInnes #endif 126304ca555eSLois Curfman McInnes } 12646d84be18SBarry Smith MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm); 126504ca555eSLois Curfman McInnes *norm = sqrt(*norm); 126604ca555eSLois Curfman McInnes } 126704ca555eSLois Curfman McInnes else if (type == NORM_1) { /* max column norm */ 126804ca555eSLois Curfman McInnes double *tmp, *tmp2; 126904ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 1270758f045eSSatish Balay tmp = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp); 1271758f045eSSatish Balay tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2); 1272cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 127304ca555eSLois Curfman McInnes *norm = 0.0; 127404ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 127504ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1276579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 127704ca555eSLois Curfman McInnes } 127804ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 127904ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1280579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 128104ca555eSLois Curfman McInnes } 12826d84be18SBarry Smith MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm); 128304ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 128404ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 128504ca555eSLois Curfman McInnes } 12860452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 128704ca555eSLois Curfman McInnes } 128804ca555eSLois Curfman McInnes else if (type == NORM_INFINITY) { /* max row norm */ 1289515d9167SLois Curfman McInnes double ntemp = 0.0; 129004ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1291dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 129204ca555eSLois Curfman McInnes sum = 0.0; 129304ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1294cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 129504ca555eSLois Curfman McInnes } 1296dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 129704ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1298cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 129904ca555eSLois Curfman McInnes } 1300515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 130104ca555eSLois Curfman McInnes } 13026d84be18SBarry Smith MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm); 130304ca555eSLois Curfman McInnes } 130404ca555eSLois Curfman McInnes else { 1305e3372554SBarry Smith SETERRQ(1,0,"No support for two norm"); 130604ca555eSLois Curfman McInnes } 130737fa93a5SLois Curfman McInnes } 1308855ac2c5SLois Curfman McInnes return 0; 1309855ac2c5SLois Curfman McInnes } 1310855ac2c5SLois Curfman McInnes 13115615d1e5SSatish Balay #undef __FUNC__ 13125615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ" 13138f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1314b7c46309SBarry Smith { 1315b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1316dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1317416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1318416022c9SBarry Smith Mat B; 1319b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 1320b7c46309SBarry Smith Scalar *array; 1321b7c46309SBarry Smith 13223638b69dSLois Curfman McInnes if (matout == PETSC_NULL && M != N) 1323e3372554SBarry Smith SETERRQ(1,0,"Square matrix only for in-place"); 1324b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,PETSC_NULL,0, 1325b4fd4287SBarry Smith PETSC_NULL,&B); CHKERRQ(ierr); 1326b7c46309SBarry Smith 1327b7c46309SBarry Smith /* copy over the A part */ 1328ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1329b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1330b7c46309SBarry Smith row = a->rstart; 1331dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1332b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1333416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1334b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1335b7c46309SBarry Smith } 1336b7c46309SBarry Smith aj = Aloc->j; 13374af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1338b7c46309SBarry Smith 1339b7c46309SBarry Smith /* copy over the B part */ 1340ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1341b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1342b7c46309SBarry Smith row = a->rstart; 13430452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1344dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1345b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1346416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1347b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1348b7c46309SBarry Smith } 13490452661fSBarry Smith PetscFree(ct); 13506d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13516d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13523638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13530de55854SLois Curfman McInnes *matout = B; 13540de55854SLois Curfman McInnes } else { 13550de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13560452661fSBarry Smith PetscFree(a->rowners); 13570de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13580de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 13590452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 13600452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13610de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1362a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13630452661fSBarry Smith PetscFree(a); 1364*f09e8eb9SSatish Balay PetscMemcpy(A,B,sizeof(struct _p_Mat)); 13650452661fSBarry Smith PetscHeaderDestroy(B); 13660de55854SLois Curfman McInnes } 1367b7c46309SBarry Smith return 0; 1368b7c46309SBarry Smith } 1369b7c46309SBarry Smith 13705615d1e5SSatish Balay #undef __FUNC__ 13715615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ" 13724b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr) 1373a008b906SSatish Balay { 13744b967eb1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 13754b967eb1SSatish Balay Mat a = aij->A, b = aij->B; 1376a008b906SSatish Balay int ierr,s1,s2,s3; 1377a008b906SSatish Balay 13784b967eb1SSatish Balay ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr); 13794b967eb1SSatish Balay if (rr) { 13804b967eb1SSatish Balay s3 = aij->n; 13814b967eb1SSatish Balay VecGetLocalSize_Fast(rr,s1); 1382e3372554SBarry Smith if (s1!=s3) SETERRQ(1,0,"right vector non-conforming local size"); 13834b967eb1SSatish Balay /* Overlap communication with computation. */ 138443a90d84SBarry Smith ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1385a008b906SSatish Balay } 13864b967eb1SSatish Balay if (ll) { 13874b967eb1SSatish Balay VecGetLocalSize_Fast(ll,s1); 1388e3372554SBarry Smith if (s1!=s2) SETERRQ(1,0,"left vector non-conforming local size"); 1389c351683dSSatish Balay ierr = (*b->ops.diagonalscale)(b,ll,0); CHKERRQ(ierr); 13904b967eb1SSatish Balay } 13914b967eb1SSatish Balay /* scale the diagonal block */ 1392c351683dSSatish Balay ierr = (*a->ops.diagonalscale)(a,ll,rr); CHKERRQ(ierr); 13934b967eb1SSatish Balay 13944b967eb1SSatish Balay if (rr) { 13954b967eb1SSatish Balay /* Do a scatter end and then right scale the off-diagonal block */ 139643a90d84SBarry Smith ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1397c351683dSSatish Balay ierr = (*b->ops.diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr); 13984b967eb1SSatish Balay } 13994b967eb1SSatish Balay 1400a008b906SSatish Balay return 0; 1401a008b906SSatish Balay } 1402a008b906SSatish Balay 1403a008b906SSatish Balay 1404682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat); 14055615d1e5SSatish Balay #undef __FUNC__ 14065eea60f9SBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" /* ADIC Ignore */ 14078f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A) 1408682d7d0cSBarry Smith { 1409682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1410682d7d0cSBarry Smith 1411682d7d0cSBarry Smith if (!a->rank) return MatPrintHelp_SeqAIJ(a->A); 1412682d7d0cSBarry Smith else return 0; 1413682d7d0cSBarry Smith } 1414682d7d0cSBarry Smith 14155615d1e5SSatish Balay #undef __FUNC__ 14165eea60f9SBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" /* ADIC Ignore */ 14178f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs) 14185a838052SSatish Balay { 14195a838052SSatish Balay *bs = 1; 14205a838052SSatish Balay return 0; 14215a838052SSatish Balay } 14225615d1e5SSatish Balay #undef __FUNC__ 14235eea60f9SBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" /* ADIC Ignore */ 14248f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A) 1425bb5a7306SBarry Smith { 1426bb5a7306SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1427bb5a7306SBarry Smith int ierr; 1428bb5a7306SBarry Smith ierr = MatSetUnfactored(a->A); CHKERRQ(ierr); 1429bb5a7306SBarry Smith return 0; 1430bb5a7306SBarry Smith } 1431bb5a7306SBarry Smith 14328f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int); 14332f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 14340a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 14350a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **); 14368a729477SBarry Smith /* -------------------------------------------------------------------*/ 14372ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ, 143839e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 143944a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 144044a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 144136ce4990SBarry Smith 0,0, 144236ce4990SBarry Smith 0,0, 144336ce4990SBarry Smith 0,0, 144444a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1445b7c46309SBarry Smith MatTranspose_MPIAIJ, 1446a66be287SLois Curfman McInnes MatGetInfo_MPIAIJ,0, 1447a008b906SSatish Balay MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ, 1448ee50ffe9SBarry Smith MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ, 14491eb62cbbSBarry Smith 0, 1450299609e3SLois Curfman McInnes MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ, 145136ce4990SBarry Smith 0,0,0,0, 1452d1710a03SLois Curfman McInnes MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ, 145336ce4990SBarry Smith 0,0, 145494a9d846SBarry Smith 0,0,MatConvertSameType_MPIAIJ,0,0, 1455b49de8d1SLois Curfman McInnes 0,0,0, 1456598137ffSSatish Balay MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0, 1457052efed2SBarry Smith MatPrintHelp_MPIAIJ, 14583b2fbd54SBarry Smith MatScale_MPIAIJ,0,0,0, 14590a198c4cSBarry Smith MatGetBlockSize_MPIAIJ,0,0,0,0, 1460bb5a7306SBarry Smith MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ}; 146136ce4990SBarry Smith 14628a729477SBarry Smith 14635615d1e5SSatish Balay #undef __FUNC__ 14645615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ" 14651987afe7SBarry Smith /*@C 1466ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 14673a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 14683a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 14693a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 14703a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 14718a729477SBarry Smith 14728a729477SBarry Smith Input Parameters: 14731eb62cbbSBarry Smith . comm - MPI communicator 14747d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 147592e8d321SLois Curfman McInnes This value should be the same as the local size used in creating the 147692e8d321SLois Curfman McInnes y vector for the matrix-vector product y = Ax. 14771a3896d6SBarry Smith . n - This value should be the same as the local size used in creating the 14781a3896d6SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 14791a3896d6SBarry Smith calculated if N is given) For square matrices n is almost always m. 14807d3e4905SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 148192e8d321SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1482ab693e5aSLois Curfman McInnes . d_nz - number of nonzeros per row in diagonal portion of local submatrix 1483ff756334SLois Curfman McInnes (same for all local rows) 14842bd5e0b2SLois Curfman McInnes . d_nzz - array containing the number of nonzeros in the various rows of the 148592e8d321SLois Curfman McInnes diagonal portion of the local submatrix (possibly different for each row) 14862bd5e0b2SLois Curfman McInnes or PETSC_NULL. You must leave room for the diagonal entry even if 14872bd5e0b2SLois Curfman McInnes it is zero. 14882bd5e0b2SLois Curfman McInnes . o_nz - number of nonzeros per row in the off-diagonal portion of local 1489ab693e5aSLois Curfman McInnes submatrix (same for all local rows). 14902bd5e0b2SLois Curfman McInnes . o_nzz - array containing the number of nonzeros in the various rows of the 14912bd5e0b2SLois Curfman McInnes off-diagonal portion of the local submatrix (possibly different for 14922bd5e0b2SLois Curfman McInnes each row) or PETSC_NULL. 14938a729477SBarry Smith 1494ff756334SLois Curfman McInnes Output Parameter: 149544cd7ae7SLois Curfman McInnes . A - the matrix 14968a729477SBarry Smith 1497ff756334SLois Curfman McInnes Notes: 1498ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1499ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 15000002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 15010002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1502ff756334SLois Curfman McInnes 1503ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1504ff756334SLois Curfman McInnes (possibly both). 1505ff756334SLois Curfman McInnes 15065511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 15075511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 15085511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 15095511cfe3SLois Curfman McInnes 15105511cfe3SLois Curfman McInnes Options Database Keys: 15115511cfe3SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 15126e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 15136e62573dSLois Curfman McInnes $ (max limit=5) 15146e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 15156e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 15166e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 15175511cfe3SLois Curfman McInnes 1518e0245417SLois Curfman McInnes Storage Information: 1519e0245417SLois Curfman McInnes For a square global matrix we define each processor's diagonal portion 1520e0245417SLois Curfman McInnes to be its local rows and the corresponding columns (a square submatrix); 1521e0245417SLois Curfman McInnes each processor's off-diagonal portion encompasses the remainder of the 1522e0245417SLois Curfman McInnes local matrix (a rectangular submatrix). 1523e0245417SLois Curfman McInnes 1524e0245417SLois Curfman McInnes The user can specify preallocated storage for the diagonal part of 15255ace5be8SLois Curfman McInnes the local submatrix with either d_nz or d_nnz (not both). Set 15265ace5be8SLois Curfman McInnes d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic 15275ace5be8SLois Curfman McInnes memory allocation. Likewise, specify preallocated storage for the 15285ace5be8SLois Curfman McInnes off-diagonal part of the local submatrix with o_nz or o_nnz (not both). 1529ff756334SLois Curfman McInnes 15305511cfe3SLois Curfman McInnes Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In 15315511cfe3SLois Curfman McInnes the figure below we depict these three local rows and all columns (0-11). 15322191d07cSBarry Smith 1533b810aeb4SBarry Smith $ 0 1 2 3 4 5 6 7 8 9 10 11 1534b810aeb4SBarry Smith $ ------------------- 15355511cfe3SLois Curfman McInnes $ row 3 | o o o d d d o o o o o o 15365511cfe3SLois Curfman McInnes $ row 4 | o o o d d d o o o o o o 15375511cfe3SLois Curfman McInnes $ row 5 | o o o d d d o o o o o o 15385511cfe3SLois Curfman McInnes $ ------------------- 1539b810aeb4SBarry Smith $ 1540b810aeb4SBarry Smith 15415511cfe3SLois Curfman McInnes Thus, any entries in the d locations are stored in the d (diagonal) 15425511cfe3SLois Curfman McInnes submatrix, and any entries in the o locations are stored in the 15435511cfe3SLois Curfman McInnes o (off-diagonal) submatrix. Note that the d and the o submatrices are 15445511cfe3SLois Curfman McInnes stored simply in the MATSEQAIJ format for compressed row storage. 15455511cfe3SLois Curfman McInnes 15465511cfe3SLois Curfman McInnes Now d_nz should indicate the number of nonzeros per row in the d matrix, 15475511cfe3SLois Curfman McInnes and o_nz should indicate the number of nonzeros per row in the o matrix. 15485511cfe3SLois Curfman McInnes In general, for PDE problems in which most nonzeros are near the diagonal, 15493d323bbdSBarry Smith one expects d_nz >> o_nz. For large problems you MUST preallocate memory 155092e8d321SLois Curfman McInnes or you will get TERRIBLE performance; see the users' manual chapter on 15516da5968aSLois Curfman McInnes matrices. 15523a511b96SLois Curfman McInnes 1553dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1554ff756334SLois Curfman McInnes 1555fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 15568a729477SBarry Smith @*/ 15571eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 155844cd7ae7SLois Curfman McInnes int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A) 15598a729477SBarry Smith { 156044cd7ae7SLois Curfman McInnes Mat B; 156144cd7ae7SLois Curfman McInnes Mat_MPIAIJ *b; 156236ce4990SBarry Smith int ierr, i,sum[2],work[2],size; 1563416022c9SBarry Smith 156444cd7ae7SLois Curfman McInnes *A = 0; 1565*f09e8eb9SSatish Balay PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATMPIAIJ,comm); 156644cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 156744cd7ae7SLois Curfman McInnes B->data = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b); 156844cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_MPIAIJ)); 156944cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 157036ce4990SBarry Smith MPI_Comm_size(comm,&size); 157136ce4990SBarry Smith if (size == 1) { 157236ce4990SBarry Smith B->ops.getrowij = MatGetRowIJ_MPIAIJ; 157336ce4990SBarry Smith B->ops.restorerowij = MatRestoreRowIJ_MPIAIJ; 157436ce4990SBarry Smith B->ops.lufactorsymbolic = MatLUFactorSymbolic_MPIAIJ; 157536ce4990SBarry Smith B->ops.lufactornumeric = MatLUFactorNumeric_MPIAIJ; 157636ce4990SBarry Smith B->ops.lufactor = MatLUFactor_MPIAIJ; 157736ce4990SBarry Smith B->ops.solve = MatSolve_MPIAIJ; 157836ce4990SBarry Smith B->ops.solveadd = MatSolveAdd_MPIAIJ; 157936ce4990SBarry Smith B->ops.solvetrans = MatSolveTrans_MPIAIJ; 158036ce4990SBarry Smith B->ops.solvetransadd = MatSolveTransAdd_MPIAIJ; 158136ce4990SBarry Smith B->ops.ilufactorsymbolic = MatILUFactorSymbolic_MPIAIJ; 158236ce4990SBarry Smith } 158344cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_MPIAIJ; 158444cd7ae7SLois Curfman McInnes B->view = MatView_MPIAIJ; 158544cd7ae7SLois Curfman McInnes B->factor = 0; 158644cd7ae7SLois Curfman McInnes B->assembled = PETSC_FALSE; 158790f02eecSBarry Smith B->mapping = 0; 1588d6dfbf8fSBarry Smith 158947794344SBarry Smith B->insertmode = NOT_SET_VALUES; 159044cd7ae7SLois Curfman McInnes MPI_Comm_rank(comm,&b->rank); 159144cd7ae7SLois Curfman McInnes MPI_Comm_size(comm,&b->size); 15921eb62cbbSBarry Smith 1593b4fd4287SBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) 1594e3372554SBarry Smith SETERRQ(1,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz"); 15951987afe7SBarry Smith 1596dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE || N == PETSC_DECIDE) { 15971eb62cbbSBarry Smith work[0] = m; work[1] = n; 1598d65a2f8fSBarry Smith MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm ); 1599dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE) M = sum[0]; 1600dbd7a890SLois Curfman McInnes if (N == PETSC_DECIDE) N = sum[1]; 16011eb62cbbSBarry Smith } 160244cd7ae7SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);} 160344cd7ae7SLois Curfman McInnes if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);} 160444cd7ae7SLois Curfman McInnes b->m = m; B->m = m; 160544cd7ae7SLois Curfman McInnes b->n = n; B->n = n; 160644cd7ae7SLois Curfman McInnes b->N = N; B->N = N; 160744cd7ae7SLois Curfman McInnes b->M = M; B->M = M; 16081eb62cbbSBarry Smith 16091eb62cbbSBarry Smith /* build local table of row and column ownerships */ 161044cd7ae7SLois Curfman McInnes b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners); 1611*f09e8eb9SSatish Balay PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1612603f58a4SSatish Balay b->cowners = b->rowners + b->size + 2; 161344cd7ae7SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm); 161444cd7ae7SLois Curfman McInnes b->rowners[0] = 0; 161544cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 161644cd7ae7SLois Curfman McInnes b->rowners[i] += b->rowners[i-1]; 16178a729477SBarry Smith } 161844cd7ae7SLois Curfman McInnes b->rstart = b->rowners[b->rank]; 161944cd7ae7SLois Curfman McInnes b->rend = b->rowners[b->rank+1]; 162044cd7ae7SLois Curfman McInnes MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm); 162144cd7ae7SLois Curfman McInnes b->cowners[0] = 0; 162244cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 162344cd7ae7SLois Curfman McInnes b->cowners[i] += b->cowners[i-1]; 16248a729477SBarry Smith } 162544cd7ae7SLois Curfman McInnes b->cstart = b->cowners[b->rank]; 162644cd7ae7SLois Curfman McInnes b->cend = b->cowners[b->rank+1]; 16278a729477SBarry Smith 16285ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 1629029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr); 163044cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->A); 16317b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 1632029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr); 163344cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->B); 16348a729477SBarry Smith 16351eb62cbbSBarry Smith /* build cache for off array entries formed */ 163644cd7ae7SLois Curfman McInnes ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr); 163790f02eecSBarry Smith b->donotstash = 0; 163844cd7ae7SLois Curfman McInnes b->colmap = 0; 163944cd7ae7SLois Curfman McInnes b->garray = 0; 164044cd7ae7SLois Curfman McInnes b->roworiented = 1; 16418a729477SBarry Smith 16421eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 164344cd7ae7SLois Curfman McInnes b->lvec = 0; 164444cd7ae7SLois Curfman McInnes b->Mvctx = 0; 16458a729477SBarry Smith 16467a0afa10SBarry Smith /* stuff for MatGetRow() */ 164744cd7ae7SLois Curfman McInnes b->rowindices = 0; 164844cd7ae7SLois Curfman McInnes b->rowvalues = 0; 164944cd7ae7SLois Curfman McInnes b->getrowactive = PETSC_FALSE; 16507a0afa10SBarry Smith 165144cd7ae7SLois Curfman McInnes *A = B; 1652d6dfbf8fSBarry Smith return 0; 1653d6dfbf8fSBarry Smith } 1654c74985f6SBarry Smith 16555615d1e5SSatish Balay #undef __FUNC__ 16565615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ" 16578f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues) 1658d6dfbf8fSBarry Smith { 1659d6dfbf8fSBarry Smith Mat mat; 1660416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 1661a1b97e82SLois Curfman McInnes int ierr, len=0, flg; 1662d6dfbf8fSBarry Smith 1663416022c9SBarry Smith *newmat = 0; 1664*f09e8eb9SSatish Balay PetscHeaderCreate(mat,_p_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm); 1665a5a9c739SBarry Smith PLogObjectCreate(mat); 16660452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1667416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 166844a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 166944a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1670d6dfbf8fSBarry Smith mat->factor = matin->factor; 1671c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1672d6dfbf8fSBarry Smith 167344cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 167444cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 167544cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 167644cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 1677d6dfbf8fSBarry Smith 1678416022c9SBarry Smith a->rstart = oldmat->rstart; 1679416022c9SBarry Smith a->rend = oldmat->rend; 1680416022c9SBarry Smith a->cstart = oldmat->cstart; 1681416022c9SBarry Smith a->cend = oldmat->cend; 168217699dbbSLois Curfman McInnes a->size = oldmat->size; 168317699dbbSLois Curfman McInnes a->rank = oldmat->rank; 168447794344SBarry Smith mat->insertmode = NOT_SET_VALUES; 1685bcd2baecSBarry Smith a->rowvalues = 0; 1686bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1687d6dfbf8fSBarry Smith 1688603f58a4SSatish Balay a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1689*f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1690603f58a4SSatish Balay a->cowners = a->rowners + a->size + 2; 1691603f58a4SSatish Balay PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int)); 1692416022c9SBarry Smith ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 16932ee70a88SLois Curfman McInnes if (oldmat->colmap) { 16940452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1695416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1696416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1697416022c9SBarry Smith } else a->colmap = 0; 16983f41c07dSBarry Smith if (oldmat->garray) { 16993f41c07dSBarry Smith len = ((Mat_SeqAIJ *) (oldmat->B->data))->n; 17003f41c07dSBarry Smith a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray); 1701464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 17023f41c07dSBarry Smith if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1703416022c9SBarry Smith } else a->garray = 0; 1704d6dfbf8fSBarry Smith 1705416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1706416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1707a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1708416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 1709416022c9SBarry Smith ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 1710416022c9SBarry Smith PLogObjectParent(mat,a->A); 1711416022c9SBarry Smith ierr = MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr); 1712416022c9SBarry Smith PLogObjectParent(mat,a->B); 17135dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 171425cdf11fSBarry Smith if (flg) { 1715682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1716682d7d0cSBarry Smith } 17178a729477SBarry Smith *newmat = mat; 17188a729477SBarry Smith return 0; 17198a729477SBarry Smith } 1720416022c9SBarry Smith 172177c4ece6SBarry Smith #include "sys.h" 1722416022c9SBarry Smith 17235615d1e5SSatish Balay #undef __FUNC__ 17245615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ" 172519bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat) 1726416022c9SBarry Smith { 1727d65a2f8fSBarry Smith Mat A; 1728416022c9SBarry Smith int i, nz, ierr, j,rstart, rend, fd; 1729d65a2f8fSBarry Smith Scalar *vals,*svals; 173019bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 1731416022c9SBarry Smith MPI_Status status; 173217699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1733d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 173419bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 1735416022c9SBarry Smith 173617699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 173717699dbbSLois Curfman McInnes if (!rank) { 173890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 173977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr); 1740e3372554SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object"); 1741416022c9SBarry Smith } 1742416022c9SBarry Smith 1743416022c9SBarry Smith MPI_Bcast(header+1,3,MPI_INT,0,comm); 1744416022c9SBarry Smith M = header[1]; N = header[2]; 1745416022c9SBarry Smith /* determine ownership of all rows */ 174617699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 17470452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1748416022c9SBarry Smith MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 1749416022c9SBarry Smith rowners[0] = 0; 175017699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1751416022c9SBarry Smith rowners[i] += rowners[i-1]; 1752416022c9SBarry Smith } 175317699dbbSLois Curfman McInnes rstart = rowners[rank]; 175417699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1755416022c9SBarry Smith 1756416022c9SBarry Smith /* distribute row lengths to all processors */ 17570452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1758416022c9SBarry Smith offlens = ourlens + (rend-rstart); 175917699dbbSLois Curfman McInnes if (!rank) { 17600452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 176177c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 17620452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 176317699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1764d65a2f8fSBarry Smith MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 17650452661fSBarry Smith PetscFree(sndcounts); 1766416022c9SBarry Smith } 1767416022c9SBarry Smith else { 1768416022c9SBarry Smith MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 1769416022c9SBarry Smith } 1770416022c9SBarry Smith 177117699dbbSLois Curfman McInnes if (!rank) { 1772416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 17730452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1774cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 177517699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1776416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1777416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1778416022c9SBarry Smith } 1779416022c9SBarry Smith } 17800452661fSBarry Smith PetscFree(rowlengths); 1781416022c9SBarry Smith 1782416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1783416022c9SBarry Smith maxnz = 0; 178417699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 17850452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1786416022c9SBarry Smith } 17870452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1788416022c9SBarry Smith 1789416022c9SBarry Smith /* read in my part of the matrix column indices */ 1790416022c9SBarry Smith nz = procsnz[0]; 17910452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 179277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr); 1793d65a2f8fSBarry Smith 1794d65a2f8fSBarry Smith /* read in every one elses and ship off */ 179517699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1796d65a2f8fSBarry Smith nz = procsnz[i]; 179777c4ece6SBarry Smith ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 1798d65a2f8fSBarry Smith MPI_Send(cols,nz,MPI_INT,i,tag,comm); 1799d65a2f8fSBarry Smith } 18000452661fSBarry Smith PetscFree(cols); 1801416022c9SBarry Smith } 1802416022c9SBarry Smith else { 1803416022c9SBarry Smith /* determine buffer space needed for message */ 1804416022c9SBarry Smith nz = 0; 1805416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1806416022c9SBarry Smith nz += ourlens[i]; 1807416022c9SBarry Smith } 18080452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1809416022c9SBarry Smith 1810416022c9SBarry Smith /* receive message of column indices*/ 1811d65a2f8fSBarry Smith MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 1812416022c9SBarry Smith MPI_Get_count(&status,MPI_INT,&maxnz); 1813e3372554SBarry Smith if (maxnz != nz) SETERRQ(1,0,"something is wrong with file"); 1814416022c9SBarry Smith } 1815416022c9SBarry Smith 1816416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 1817cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 1818416022c9SBarry Smith jj = 0; 1819416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1820416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 1821d65a2f8fSBarry Smith if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1822416022c9SBarry Smith jj++; 1823416022c9SBarry Smith } 1824416022c9SBarry Smith } 1825d65a2f8fSBarry Smith 1826d65a2f8fSBarry Smith /* create our matrix */ 1827416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1828416022c9SBarry Smith ourlens[i] -= offlens[i]; 1829416022c9SBarry Smith } 1830d65a2f8fSBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 1831d65a2f8fSBarry Smith A = *newmat; 18326d4a8577SBarry Smith MatSetOption(A,MAT_COLUMNS_SORTED); 1833d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1834d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 1835d65a2f8fSBarry Smith } 1836416022c9SBarry Smith 183717699dbbSLois Curfman McInnes if (!rank) { 18380452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1839416022c9SBarry Smith 1840416022c9SBarry Smith /* read in my part of the matrix numerical values */ 1841416022c9SBarry Smith nz = procsnz[0]; 184277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 1843d65a2f8fSBarry Smith 1844d65a2f8fSBarry Smith /* insert into matrix */ 1845d65a2f8fSBarry Smith jj = rstart; 1846d65a2f8fSBarry Smith smycols = mycols; 1847d65a2f8fSBarry Smith svals = vals; 1848d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1849d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1850d65a2f8fSBarry Smith smycols += ourlens[i]; 1851d65a2f8fSBarry Smith svals += ourlens[i]; 1852d65a2f8fSBarry Smith jj++; 1853416022c9SBarry Smith } 1854416022c9SBarry Smith 1855d65a2f8fSBarry Smith /* read in other processors and ship out */ 185617699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1857416022c9SBarry Smith nz = procsnz[i]; 185877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 1859d65a2f8fSBarry Smith MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 1860416022c9SBarry Smith } 18610452661fSBarry Smith PetscFree(procsnz); 1862416022c9SBarry Smith } 1863d65a2f8fSBarry Smith else { 1864d65a2f8fSBarry Smith /* receive numeric values */ 18650452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1866416022c9SBarry Smith 1867d65a2f8fSBarry Smith /* receive message of values*/ 1868d65a2f8fSBarry Smith MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 1869d65a2f8fSBarry Smith MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 1870e3372554SBarry Smith if (maxnz != nz) SETERRQ(1,0,"something is wrong with file"); 1871d65a2f8fSBarry Smith 1872d65a2f8fSBarry Smith /* insert into matrix */ 1873d65a2f8fSBarry Smith jj = rstart; 1874d65a2f8fSBarry Smith smycols = mycols; 1875d65a2f8fSBarry Smith svals = vals; 1876d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1877d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1878d65a2f8fSBarry Smith smycols += ourlens[i]; 1879d65a2f8fSBarry Smith svals += ourlens[i]; 1880d65a2f8fSBarry Smith jj++; 1881d65a2f8fSBarry Smith } 1882d65a2f8fSBarry Smith } 18830452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1884d65a2f8fSBarry Smith 18856d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 18866d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1887416022c9SBarry Smith return 0; 1888416022c9SBarry Smith } 1889