1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*5a655dc6SBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.289 1999/03/18 00:35:56 balay Exp bsmith $"; 3cb512458SBarry Smith #endif 48a729477SBarry Smith 570f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 6f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 7d9942c19SSatish Balay #include "src/inline/spops.h" 88a729477SBarry Smith 9bc5ccf88SSatish Balay extern int MatSetUpMultiply_MPIAIJ(Mat); 10bc5ccf88SSatish Balay extern int DisAssemble_MPIAIJ(Mat); 11bc5ccf88SSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode); 12bc5ccf88SSatish Balay extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 13bc5ccf88SSatish Balay extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 14bc5ccf88SSatish Balay extern int MatPrintHelp_SeqAIJ(Mat); 15bc5ccf88SSatish Balay 169e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column 179e25ed09SBarry Smith number to the local number in the off-diagonal part of the local 189e25ed09SBarry Smith storage of the matrix. This is done in a non scable way since the 199e25ed09SBarry Smith length of colmap equals the global matrix length. 209e25ed09SBarry Smith */ 215615d1e5SSatish Balay #undef __FUNC__ 22d4bb536fSBarry Smith #define __FUNC__ "CreateColmap_MPIAIJ_Private" 230a198c4cSBarry Smith int CreateColmap_MPIAIJ_Private(Mat mat) 249e25ed09SBarry Smith { 2544a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 26ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data; 276dee3f7eSBarry Smith int n = B->n,i; 286dee3f7eSBarry Smith #if defined (USE_CTABLE) 296dee3f7eSBarry Smith int ierr; 306dee3f7eSBarry Smith #endif 31dbb450caSBarry Smith 323a40ed3dSBarry Smith PetscFunctionBegin; 33b1fc9764SSatish Balay #if defined (USE_CTABLE) 34fa46199cSSatish Balay ierr = TableCreate(aij->n/5,&aij->colmap); CHKERRQ(ierr); 35b1fc9764SSatish Balay for ( i=0; i<n; i++ ){ 36b1fc9764SSatish Balay ierr = TableAdd(aij->colmap,aij->garray[i]+1,i+1);CHKERRQ(ierr); 37b1fc9764SSatish Balay } 38b1fc9764SSatish Balay #else 39758f045eSSatish Balay aij->colmap = (int *) PetscMalloc((aij->N+1)*sizeof(int));CHKPTRQ(aij->colmap); 40464493b3SBarry Smith PLogObjectMemory(mat,aij->N*sizeof(int)); 41cddf8d76SBarry Smith PetscMemzero(aij->colmap,aij->N*sizeof(int)); 42905e6a2fSBarry Smith for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1; 43b1fc9764SSatish Balay #endif 443a40ed3dSBarry Smith PetscFunctionReturn(0); 459e25ed09SBarry Smith } 469e25ed09SBarry Smith 470520107fSSatish Balay #define CHUNKSIZE 15 4830770e4dSSatish Balay #define MatSetValues_SeqAIJ_A_Private(row,col,value,addv) \ 490520107fSSatish Balay { \ 500520107fSSatish Balay \ 510520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 5230770e4dSSatish Balay rmax = aimax[row]; nrow = ailen[row]; \ 53f5e9677aSSatish Balay col1 = col - shift; \ 54f5e9677aSSatish Balay \ 55ba4e3ef2SSatish Balay low = 0; high = nrow; \ 56ba4e3ef2SSatish Balay while (high-low > 5) { \ 57ba4e3ef2SSatish Balay t = (low+high)/2; \ 58ba4e3ef2SSatish Balay if (rp[t] > col) high = t; \ 59ba4e3ef2SSatish Balay else low = t; \ 60ba4e3ef2SSatish Balay } \ 610520107fSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 62f5e9677aSSatish Balay if (rp[_i] > col1) break; \ 63f5e9677aSSatish Balay if (rp[_i] == col1) { \ 640520107fSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 650520107fSSatish Balay else ap[_i] = value; \ 6630770e4dSSatish Balay goto a_noinsert; \ 670520107fSSatish Balay } \ 680520107fSSatish Balay } \ 6989280ab3SLois Curfman McInnes if (nonew == 1) goto a_noinsert; \ 70a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero into matrix"); \ 710520107fSSatish Balay if (nrow >= rmax) { \ 720520107fSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 730520107fSSatish Balay int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \ 740520107fSSatish Balay Scalar *new_a; \ 750520107fSSatish Balay \ 76a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); \ 7789280ab3SLois Curfman McInnes \ 780520107fSSatish Balay /* malloc new storage space */ \ 790520107fSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \ 800520107fSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 810520107fSSatish Balay new_j = (int *) (new_a + new_nz); \ 820520107fSSatish Balay new_i = new_j + new_nz; \ 830520107fSSatish Balay \ 840520107fSSatish Balay /* copy over old data into new slots */ \ 850520107fSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \ 860520107fSSatish Balay for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \ 870520107fSSatish Balay PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \ 880520107fSSatish Balay len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \ 890520107fSSatish Balay PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \ 900520107fSSatish Balay len*sizeof(int)); \ 910520107fSSatish Balay PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \ 920520107fSSatish Balay PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \ 930520107fSSatish Balay len*sizeof(Scalar)); \ 940520107fSSatish Balay /* free up old matrix storage */ \ 95f5e9677aSSatish Balay \ 960520107fSSatish Balay PetscFree(a->a); \ 970520107fSSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \ 980520107fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; \ 990520107fSSatish Balay a->singlemalloc = 1; \ 1000520107fSSatish Balay \ 1010520107fSSatish Balay rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; \ 10230770e4dSSatish Balay rmax = aimax[row] = aimax[row] + CHUNKSIZE; \ 1030520107fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 1040520107fSSatish Balay a->maxnz += CHUNKSIZE; \ 1050520107fSSatish Balay a->reallocs++; \ 1060520107fSSatish Balay } \ 1070520107fSSatish Balay N = nrow++ - 1; a->nz++; \ 1080520107fSSatish Balay /* shift up all the later entries in this row */ \ 1090520107fSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 1100520107fSSatish Balay rp[ii+1] = rp[ii]; \ 1110520107fSSatish Balay ap[ii+1] = ap[ii]; \ 1120520107fSSatish Balay } \ 113f5e9677aSSatish Balay rp[_i] = col1; \ 1140520107fSSatish Balay ap[_i] = value; \ 11530770e4dSSatish Balay a_noinsert: ; \ 1160520107fSSatish Balay ailen[row] = nrow; \ 1170520107fSSatish Balay } 1180a198c4cSBarry Smith 11930770e4dSSatish Balay #define MatSetValues_SeqAIJ_B_Private(row,col,value,addv) \ 12030770e4dSSatish Balay { \ 12130770e4dSSatish Balay \ 12230770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 12330770e4dSSatish Balay rmax = bimax[row]; nrow = bilen[row]; \ 12430770e4dSSatish Balay col1 = col - shift; \ 12530770e4dSSatish Balay \ 126ba4e3ef2SSatish Balay low = 0; high = nrow; \ 127ba4e3ef2SSatish Balay while (high-low > 5) { \ 128ba4e3ef2SSatish Balay t = (low+high)/2; \ 129ba4e3ef2SSatish Balay if (rp[t] > col) high = t; \ 130ba4e3ef2SSatish Balay else low = t; \ 131ba4e3ef2SSatish Balay } \ 13230770e4dSSatish Balay for ( _i=0; _i<nrow; _i++ ) { \ 13330770e4dSSatish Balay if (rp[_i] > col1) break; \ 13430770e4dSSatish Balay if (rp[_i] == col1) { \ 13530770e4dSSatish Balay if (addv == ADD_VALUES) ap[_i] += value; \ 13630770e4dSSatish Balay else ap[_i] = value; \ 13730770e4dSSatish Balay goto b_noinsert; \ 13830770e4dSSatish Balay } \ 13930770e4dSSatish Balay } \ 14089280ab3SLois Curfman McInnes if (nonew == 1) goto b_noinsert; \ 141a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero into matrix"); \ 14230770e4dSSatish Balay if (nrow >= rmax) { \ 14330770e4dSSatish Balay /* there is no extra room in row, therefore enlarge */ \ 14474c639caSSatish Balay int new_nz = bi[b->m] + CHUNKSIZE,len,*new_i,*new_j; \ 14530770e4dSSatish Balay Scalar *new_a; \ 14630770e4dSSatish Balay \ 147a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); \ 14889280ab3SLois Curfman McInnes \ 14930770e4dSSatish Balay /* malloc new storage space */ \ 15074c639caSSatish Balay len = new_nz*(sizeof(int)+sizeof(Scalar))+(b->m+1)*sizeof(int); \ 15130770e4dSSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \ 15230770e4dSSatish Balay new_j = (int *) (new_a + new_nz); \ 15330770e4dSSatish Balay new_i = new_j + new_nz; \ 15430770e4dSSatish Balay \ 15530770e4dSSatish Balay /* copy over old data into new slots */ \ 15630770e4dSSatish Balay for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = bi[ii];} \ 15774c639caSSatish Balay for ( ii=row+1; ii<b->m+1; ii++ ) {new_i[ii] = bi[ii]+CHUNKSIZE;} \ 15830770e4dSSatish Balay PetscMemcpy(new_j,bj,(bi[row]+nrow+shift)*sizeof(int)); \ 15930770e4dSSatish Balay len = (new_nz - CHUNKSIZE - bi[row] - nrow - shift); \ 16030770e4dSSatish Balay PetscMemcpy(new_j+bi[row]+shift+nrow+CHUNKSIZE,bj+bi[row]+shift+nrow, \ 16130770e4dSSatish Balay len*sizeof(int)); \ 16230770e4dSSatish Balay PetscMemcpy(new_a,ba,(bi[row]+nrow+shift)*sizeof(Scalar)); \ 16330770e4dSSatish Balay PetscMemcpy(new_a+bi[row]+shift+nrow+CHUNKSIZE,ba+bi[row]+shift+nrow, \ 16430770e4dSSatish Balay len*sizeof(Scalar)); \ 16530770e4dSSatish Balay /* free up old matrix storage */ \ 16630770e4dSSatish Balay \ 16774c639caSSatish Balay PetscFree(b->a); \ 16874c639caSSatish Balay if (!b->singlemalloc) {PetscFree(b->i);PetscFree(b->j);} \ 16974c639caSSatish Balay ba = b->a = new_a; bi = b->i = new_i; bj = b->j = new_j; \ 17074c639caSSatish Balay b->singlemalloc = 1; \ 17130770e4dSSatish Balay \ 17230770e4dSSatish Balay rp = bj + bi[row] + shift; ap = ba + bi[row] + shift; \ 17330770e4dSSatish Balay rmax = bimax[row] = bimax[row] + CHUNKSIZE; \ 17474c639caSSatish Balay PLogObjectMemory(B,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \ 17574c639caSSatish Balay b->maxnz += CHUNKSIZE; \ 17674c639caSSatish Balay b->reallocs++; \ 17730770e4dSSatish Balay } \ 17874c639caSSatish Balay N = nrow++ - 1; b->nz++; \ 17930770e4dSSatish Balay /* shift up all the later entries in this row */ \ 18030770e4dSSatish Balay for ( ii=N; ii>=_i; ii-- ) { \ 18130770e4dSSatish Balay rp[ii+1] = rp[ii]; \ 18230770e4dSSatish Balay ap[ii+1] = ap[ii]; \ 18330770e4dSSatish Balay } \ 18430770e4dSSatish Balay rp[_i] = col1; \ 18530770e4dSSatish Balay ap[_i] = value; \ 18630770e4dSSatish Balay b_noinsert: ; \ 18730770e4dSSatish Balay bilen[row] = nrow; \ 18830770e4dSSatish Balay } 18930770e4dSSatish Balay 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 211ba4e3ef2SSatish Balay int *rp,ii,nrow,_i,rmax, N, col1,low,high,t; 21230770e4dSSatish Balay int nonew = a->nonew,shift = a->indexshift; 21330770e4dSSatish Balay Scalar *ap; 2144ee7247eSSatish Balay 2153a40ed3dSBarry Smith PetscFunctionBegin; 2168a729477SBarry Smith for ( i=0; i<m; i++ ) { 2175ef9f2a5SBarry Smith if (im[i] < 0) continue; 2183a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 219a8c6a408SBarry Smith if (im[i] >= aij->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 2200a198c4cSBarry Smith #endif 2214b0e389bSBarry Smith if (im[i] >= rstart && im[i] < rend) { 2224b0e389bSBarry Smith row = im[i] - rstart; 2231eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 2244b0e389bSBarry Smith if (in[j] >= cstart && in[j] < cend){ 2254b0e389bSBarry Smith col = in[j] - cstart; 2264b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 22730770e4dSSatish Balay MatSetValues_SeqAIJ_A_Private(row,col,value,addv); 2280520107fSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2291eb62cbbSBarry Smith } 2305ef9f2a5SBarry Smith else if (in[j] < 0) continue; 2313a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 232a8c6a408SBarry Smith else if (in[j] >= aij->N) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");} 2330a198c4cSBarry Smith #endif 2341eb62cbbSBarry Smith else { 235227d817aSBarry Smith if (mat->was_assembled) { 236905e6a2fSBarry Smith if (!aij->colmap) { 237905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 238905e6a2fSBarry Smith } 239b1fc9764SSatish Balay #if defined (USE_CTABLE) 240fa46199cSSatish Balay ierr = TableFind(aij->colmap,in[j]+1,&col); CHKERRQ(ierr); 241fa46199cSSatish Balay col--; 242b1fc9764SSatish Balay #else 243905e6a2fSBarry Smith col = aij->colmap[in[j]] - 1; 244b1fc9764SSatish Balay #endif 245ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 2462493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 2474b0e389bSBarry Smith col = in[j]; 2489bf004c3SSatish Balay /* Reinitialize the variables required by MatSetValues_SeqAIJ_B_Private() */ 249f9508a3cSSatish Balay B = aij->B; 250f9508a3cSSatish Balay b = (Mat_SeqAIJ *) B->data; 251f9508a3cSSatish Balay bimax = b->imax; bi = b->i; bilen = b->ilen; bj = b->j; 252f9508a3cSSatish Balay ba = b->a; 253d6dfbf8fSBarry Smith } 254c48de900SBarry Smith } else col = in[j]; 2554b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 25630770e4dSSatish Balay MatSetValues_SeqAIJ_B_Private(row,col,value,addv); 25730770e4dSSatish Balay /* ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */ 2581eb62cbbSBarry Smith } 2591eb62cbbSBarry Smith } 2605ef9f2a5SBarry Smith } else { 26190f02eecSBarry Smith if (!aij->donotstash) { 262d36fbae8SSatish Balay if (roworiented) { 2638798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n);CHKERRQ(ierr); 264d36fbae8SSatish Balay } else { 2658798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m);CHKERRQ(ierr); 2664b0e389bSBarry Smith } 2671eb62cbbSBarry Smith } 2688a729477SBarry Smith } 26990f02eecSBarry Smith } 2703a40ed3dSBarry Smith PetscFunctionReturn(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 2813a40ed3dSBarry Smith PetscFunctionBegin; 282b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 283a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 284a8c6a408SBarry Smith if (idxm[i] >= aij->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 285b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 286b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 287b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 288a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 289a8c6a408SBarry Smith if (idxn[j] >= aij->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 290b49de8d1SLois Curfman McInnes if (idxn[j] >= cstart && idxn[j] < cend){ 291b49de8d1SLois Curfman McInnes col = idxn[j] - cstart; 292b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 293fa852ad4SSatish Balay } else { 294905e6a2fSBarry Smith if (!aij->colmap) { 295905e6a2fSBarry Smith ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 296905e6a2fSBarry Smith } 297b1fc9764SSatish Balay #if defined (USE_CTABLE) 298fa46199cSSatish Balay ierr = TableFind(aij->colmap,idxn[j]+1,&col); CHKERRQ(ierr); 299fa46199cSSatish Balay col --; 300b1fc9764SSatish Balay #else 301905e6a2fSBarry Smith col = aij->colmap[idxn[j]] - 1; 302b1fc9764SSatish Balay #endif 303e60e1c95SSatish Balay if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0; 304d9d09a02SSatish Balay else { 305b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 306b49de8d1SLois Curfman McInnes } 307b49de8d1SLois Curfman McInnes } 308b49de8d1SLois Curfman McInnes } 309a8c6a408SBarry Smith } else { 310a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 311b49de8d1SLois Curfman McInnes } 312b49de8d1SLois Curfman McInnes } 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 314b49de8d1SLois Curfman McInnes } 315bc5ccf88SSatish Balay 316bc5ccf88SSatish Balay #undef __FUNC__ 317bc5ccf88SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ" 318bc5ccf88SSatish Balay int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode) 319bc5ccf88SSatish Balay { 320bc5ccf88SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 321d36fbae8SSatish Balay int ierr,nstash,reallocs; 322bc5ccf88SSatish Balay InsertMode addv; 323bc5ccf88SSatish Balay 324bc5ccf88SSatish Balay PetscFunctionBegin; 325bc5ccf88SSatish Balay if (aij->donotstash) { 326bc5ccf88SSatish Balay PetscFunctionReturn(0); 327bc5ccf88SSatish Balay } 328bc5ccf88SSatish Balay 329bc5ccf88SSatish Balay /* make sure all processors are either in INSERTMODE or ADDMODE */ 330bc5ccf88SSatish Balay ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,mat->comm);CHKERRQ(ierr); 331bc5ccf88SSatish Balay if (addv == (ADD_VALUES|INSERT_VALUES)) { 332bc5ccf88SSatish Balay SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Some processors inserted others added"); 333bc5ccf88SSatish Balay } 334bc5ccf88SSatish Balay mat->insertmode = addv; /* in case this processor had no cache */ 335bc5ccf88SSatish Balay 3368798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,aij->rowners); CHKERRQ(ierr); 3378798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs); CHKERRQ(ierr); 338*5a655dc6SBarry Smith PLogInfo(aij->A,"MatAssemblyBegin_MPIAIJ:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs); 339bc5ccf88SSatish Balay PetscFunctionReturn(0); 340bc5ccf88SSatish Balay } 341bc5ccf88SSatish Balay 342bc5ccf88SSatish Balay 343bc5ccf88SSatish Balay #undef __FUNC__ 344bc5ccf88SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ" 345bc5ccf88SSatish Balay int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode) 346bc5ccf88SSatish Balay { 347bc5ccf88SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 348a2d1c673SSatish Balay int i,j,rstart,ncols,n,ierr,flg; 349bc5ccf88SSatish Balay int *row,*col,other_disassembled; 350bc5ccf88SSatish Balay Scalar *val; 351bc5ccf88SSatish Balay InsertMode addv = mat->insertmode; 352bc5ccf88SSatish Balay 353bc5ccf88SSatish Balay PetscFunctionBegin; 354bc5ccf88SSatish Balay if (!aij->donotstash) { 355a2d1c673SSatish Balay while (1) { 3568798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg); CHKERRQ(ierr); 357a2d1c673SSatish Balay if (!flg) break; 358a2d1c673SSatish Balay 359bc5ccf88SSatish Balay for ( i=0; i<n; ) { 360bc5ccf88SSatish Balay /* Now identify the consecutive vals belonging to the same row */ 361bc5ccf88SSatish Balay for ( j=i,rstart=row[j]; j<n; j++ ) { if (row[j] != rstart) break; } 362bc5ccf88SSatish Balay if (j < n) ncols = j-i; 363bc5ccf88SSatish Balay else ncols = n-i; 364bc5ccf88SSatish Balay /* Now assemble all these values with a single function call */ 365bc5ccf88SSatish Balay ierr = MatSetValues_MPIAIJ(mat,1,row+i,ncols,col+i,val+i,addv); CHKERRQ(ierr); 366bc5ccf88SSatish Balay i = j; 367bc5ccf88SSatish Balay } 368bc5ccf88SSatish Balay } 3698798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash); CHKERRQ(ierr); 370bc5ccf88SSatish Balay } 371bc5ccf88SSatish Balay 372bc5ccf88SSatish Balay ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr); 373bc5ccf88SSatish Balay ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr); 374bc5ccf88SSatish Balay 375bc5ccf88SSatish Balay /* determine if any processor has disassembled, if so we must 376bc5ccf88SSatish Balay also disassemble ourselfs, in order that we may reassemble. */ 377bc5ccf88SSatish Balay /* 378bc5ccf88SSatish Balay if nonzero structure of submatrix B cannot change then we know that 379bc5ccf88SSatish Balay no processor disassembled thus we can skip this stuff 380bc5ccf88SSatish Balay */ 381bc5ccf88SSatish Balay if (!((Mat_SeqAIJ*) aij->B->data)->nonew) { 382bc5ccf88SSatish Balay ierr = MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);CHKERRQ(ierr); 383bc5ccf88SSatish Balay if (mat->was_assembled && !other_disassembled) { 384bc5ccf88SSatish Balay ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 385bc5ccf88SSatish Balay } 386bc5ccf88SSatish Balay } 387bc5ccf88SSatish Balay 388bc5ccf88SSatish Balay if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 389bc5ccf88SSatish Balay ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr); 390bc5ccf88SSatish Balay } 391bc5ccf88SSatish Balay ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr); 392bc5ccf88SSatish Balay ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr); 393bc5ccf88SSatish Balay 394bc5ccf88SSatish Balay if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;} 395bc5ccf88SSatish Balay PetscFunctionReturn(0); 396bc5ccf88SSatish Balay } 397bc5ccf88SSatish Balay 3985615d1e5SSatish Balay #undef __FUNC__ 3995615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ" 4008f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A) 4011eb62cbbSBarry Smith { 40244a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 403dbd7a890SLois Curfman McInnes int ierr; 4043a40ed3dSBarry Smith 4053a40ed3dSBarry Smith PetscFunctionBegin; 40678b31e54SBarry Smith ierr = MatZeroEntries(l->A); CHKERRQ(ierr); 40778b31e54SBarry Smith ierr = MatZeroEntries(l->B); CHKERRQ(ierr); 4083a40ed3dSBarry Smith PetscFunctionReturn(0); 4091eb62cbbSBarry Smith } 4101eb62cbbSBarry Smith 4115615d1e5SSatish Balay #undef __FUNC__ 4125615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ" 4138f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 4141eb62cbbSBarry Smith { 41544a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 41617699dbbSLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 4176a5c57faSSatish Balay int *procs,*nprocs,j,found,idx,nsends,*work,row; 41817699dbbSLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 4195392566eSBarry Smith int *rvalues,tag = A->tag,count,base,slen,n,*source; 4206a5c57faSSatish Balay int *lens,imdex,*lrows,*values,rstart=l->rstart; 421d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 4221eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 4231eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 4241eb62cbbSBarry Smith IS istmp; 4251eb62cbbSBarry Smith 4263a40ed3dSBarry Smith PetscFunctionBegin; 42777c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 42878b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 4291eb62cbbSBarry Smith 4301eb62cbbSBarry Smith /* first count number of contributors to each processor */ 4310452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 432cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 4330452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 4341eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 4351eb62cbbSBarry Smith idx = rows[i]; 4361eb62cbbSBarry Smith found = 0; 43717699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 4381eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 4391eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 4401eb62cbbSBarry Smith } 4411eb62cbbSBarry Smith } 442a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 4431eb62cbbSBarry Smith } 44417699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 4451eb62cbbSBarry Smith 4461eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 4470452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 448ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 44917699dbbSLois Curfman McInnes nrecvs = work[rank]; 450ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 45117699dbbSLois Curfman McInnes nmax = work[rank]; 4520452661fSBarry Smith PetscFree(work); 4531eb62cbbSBarry Smith 4541eb62cbbSBarry Smith /* post receives: */ 4553a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 456ca161407SBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 4571eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 458ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 4591eb62cbbSBarry Smith } 4601eb62cbbSBarry Smith 4611eb62cbbSBarry Smith /* do sends: 4621eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 4631eb62cbbSBarry Smith the ith processor 4641eb62cbbSBarry Smith */ 4650452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 4663a40ed3dSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 4670452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 4681eb62cbbSBarry Smith starts[0] = 0; 46917699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 4701eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 4711eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 4721eb62cbbSBarry Smith } 4731eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 4741eb62cbbSBarry Smith 4751eb62cbbSBarry Smith starts[0] = 0; 47617699dbbSLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 4771eb62cbbSBarry Smith count = 0; 47817699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 4791eb62cbbSBarry Smith if (procs[i]) { 480ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 4811eb62cbbSBarry Smith } 4821eb62cbbSBarry Smith } 4830452661fSBarry Smith PetscFree(starts); 4841eb62cbbSBarry Smith 48517699dbbSLois Curfman McInnes base = owners[rank]; 4861eb62cbbSBarry Smith 4871eb62cbbSBarry Smith /* wait on receives */ 4880452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 4891eb62cbbSBarry Smith source = lens + nrecvs; 4901eb62cbbSBarry Smith count = nrecvs; slen = 0; 4911eb62cbbSBarry Smith while (count) { 492ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4931eb62cbbSBarry Smith /* unpack receives into our local space */ 494ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 495d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 496d6dfbf8fSBarry Smith lens[imdex] = n; 4971eb62cbbSBarry Smith slen += n; 4981eb62cbbSBarry Smith count--; 4991eb62cbbSBarry Smith } 5000452661fSBarry Smith PetscFree(recv_waits); 5011eb62cbbSBarry Smith 5021eb62cbbSBarry Smith /* move the data into the send scatter */ 5030452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 5041eb62cbbSBarry Smith count = 0; 5051eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 5061eb62cbbSBarry Smith values = rvalues + i*nmax; 5071eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 5081eb62cbbSBarry Smith lrows[count++] = values[j] - base; 5091eb62cbbSBarry Smith } 5101eb62cbbSBarry Smith } 5110452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 5120452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 5131eb62cbbSBarry Smith 5141eb62cbbSBarry Smith /* actually zap the local rows */ 515029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 516464493b3SBarry Smith PLogObjectParent(A,istmp); 5176a5c57faSSatish Balay 5186eb55b6aSBarry Smith /* 5196eb55b6aSBarry Smith Zero the required rows. If the "diagonal block" of the matrix 5206eb55b6aSBarry Smith is square and the user wishes to set the diagonal we use seperate 5216eb55b6aSBarry Smith code so that MatSetValues() is not called for each diagonal allocating 5226eb55b6aSBarry Smith new memory, thus calling lots of mallocs and slowing things down. 5236eb55b6aSBarry Smith 5246eb55b6aSBarry Smith Contributed by: Mathew Knepley 5256eb55b6aSBarry Smith */ 526e2d53e46SBarry Smith /* must zero l->B before l->A because the (diag) case below may put values into l->B*/ 527e2d53e46SBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr); 5286eb55b6aSBarry Smith if (diag && (l->A->M == l->A->N)) { 5296eb55b6aSBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 530e2d53e46SBarry Smith } else if (diag) { 531e2d53e46SBarry Smith ierr = MatZeroRows(l->A,istmp,0); CHKERRQ(ierr); 532fa46199cSSatish Balay if (((Mat_SeqAIJ*)l->A->data)->nonew) { 533fa46199cSSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MatZeroRows() on rectangular matrices cannot be used with the Mat options\n\ 534fa46199cSSatish Balay MAT_NO_NEW_NONZERO_LOCATIONS,MAT_NEW_NONZERO_LOCATION_ERR,MAT_NEW_NONZERO_ALLOCATION_ERR"); 5356525c446SSatish Balay } 536e2d53e46SBarry Smith for ( i = 0; i < slen; i++ ) { 537e2d53e46SBarry Smith row = lrows[i] + rstart; 538e2d53e46SBarry Smith ierr = MatSetValues(A,1,&row,1,&row,diag,INSERT_VALUES);CHKERRQ(ierr); 539e2d53e46SBarry Smith } 540e2d53e46SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 541e2d53e46SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5426eb55b6aSBarry Smith } else { 5436a5c57faSSatish Balay ierr = MatZeroRows(l->A,istmp,0); CHKERRQ(ierr); 5446eb55b6aSBarry Smith } 54578b31e54SBarry Smith ierr = ISDestroy(istmp); CHKERRQ(ierr); 5466a5c57faSSatish Balay PetscFree(lrows); 54772dacd9aSBarry Smith 5481eb62cbbSBarry Smith /* wait on sends */ 5491eb62cbbSBarry Smith if (nsends) { 550ca161407SBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 551ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 5520452661fSBarry Smith PetscFree(send_status); 5531eb62cbbSBarry Smith } 5540452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 5551eb62cbbSBarry Smith 5563a40ed3dSBarry Smith PetscFunctionReturn(0); 5571eb62cbbSBarry Smith } 5581eb62cbbSBarry Smith 5595615d1e5SSatish Balay #undef __FUNC__ 5605615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ" 5618f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy) 5621eb62cbbSBarry Smith { 563416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 564fbd6ef76SBarry Smith int ierr,nt; 565416022c9SBarry Smith 5663a40ed3dSBarry Smith PetscFunctionBegin; 567a2ce50c7SBarry Smith ierr = VecGetLocalSize(xx,&nt); CHKERRQ(ierr); 568fbd6ef76SBarry Smith if (nt != a->n) { 569a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Incompatible partition of A and xx"); 570fbd6ef76SBarry Smith } 57143a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 572f830108cSBarry Smith ierr = (*a->A->ops->mult)(a->A,xx,yy); CHKERRQ(ierr); 57343a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr); 574f830108cSBarry Smith ierr = (*a->B->ops->multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr); 5753a40ed3dSBarry Smith PetscFunctionReturn(0); 5761eb62cbbSBarry Smith } 5771eb62cbbSBarry Smith 5785615d1e5SSatish Balay #undef __FUNC__ 5795615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ" 5808f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 581da3a660dSBarry Smith { 582416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 583da3a660dSBarry Smith int ierr; 5843a40ed3dSBarry Smith 5853a40ed3dSBarry Smith PetscFunctionBegin; 58643a90d84SBarry Smith ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 587f830108cSBarry Smith ierr = (*a->A->ops->multadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 58843a90d84SBarry Smith ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 589f830108cSBarry Smith ierr = (*a->B->ops->multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr); 5903a40ed3dSBarry Smith PetscFunctionReturn(0); 591da3a660dSBarry Smith } 592da3a660dSBarry Smith 5935615d1e5SSatish Balay #undef __FUNC__ 5945615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ" 5958f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy) 596da3a660dSBarry Smith { 597416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 598da3a660dSBarry Smith int ierr; 599da3a660dSBarry Smith 6003a40ed3dSBarry Smith PetscFunctionBegin; 601da3a660dSBarry Smith /* do nondiagonal part */ 602f830108cSBarry Smith ierr = (*a->B->ops->multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 603da3a660dSBarry Smith /* send it on its way */ 604537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 605da3a660dSBarry Smith /* do local part */ 606f830108cSBarry Smith ierr = (*a->A->ops->multtrans)(a->A,xx,yy); CHKERRQ(ierr); 607da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 608da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 609da3a660dSBarry Smith /* but is not perhaps always true. */ 610537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 612da3a660dSBarry Smith } 613da3a660dSBarry Smith 6145615d1e5SSatish Balay #undef __FUNC__ 6155615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ" 6168f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 617da3a660dSBarry Smith { 618416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 619da3a660dSBarry Smith int ierr; 620da3a660dSBarry Smith 6213a40ed3dSBarry Smith PetscFunctionBegin; 622da3a660dSBarry Smith /* do nondiagonal part */ 623f830108cSBarry Smith ierr = (*a->B->ops->multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr); 624da3a660dSBarry Smith /* send it on its way */ 625537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 626da3a660dSBarry Smith /* do local part */ 627f830108cSBarry Smith ierr = (*a->A->ops->multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr); 628da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 629da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 630da3a660dSBarry Smith /* but is not perhaps always true. */ 6310a198c4cSBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 6323a40ed3dSBarry Smith PetscFunctionReturn(0); 633da3a660dSBarry Smith } 634da3a660dSBarry Smith 6351eb62cbbSBarry Smith /* 6361eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 6371eb62cbbSBarry Smith diagonal block 6381eb62cbbSBarry Smith */ 6395615d1e5SSatish Balay #undef __FUNC__ 6405615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ" 6418f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v) 6421eb62cbbSBarry Smith { 6433a40ed3dSBarry Smith int ierr; 644416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 6453a40ed3dSBarry Smith 6463a40ed3dSBarry Smith PetscFunctionBegin; 647a8c6a408SBarry Smith if (a->M != a->N) SETERRQ(PETSC_ERR_SUP,0,"Supports only square matrix where A->A is diag block"); 6485baf8537SBarry Smith if (a->rstart != a->cstart || a->rend != a->cend) { 649a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"row partition must equal col partition"); 6503a40ed3dSBarry Smith } 6513a40ed3dSBarry Smith ierr = MatGetDiagonal(a->A,v);CHKERRQ(ierr); 6523a40ed3dSBarry Smith PetscFunctionReturn(0); 6531eb62cbbSBarry Smith } 6541eb62cbbSBarry Smith 6555615d1e5SSatish Balay #undef __FUNC__ 6565615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ" 6578f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A) 658052efed2SBarry Smith { 659052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 660052efed2SBarry Smith int ierr; 6613a40ed3dSBarry Smith 6623a40ed3dSBarry Smith PetscFunctionBegin; 663052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 664052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 6653a40ed3dSBarry Smith PetscFunctionReturn(0); 666052efed2SBarry Smith } 667052efed2SBarry Smith 6685615d1e5SSatish Balay #undef __FUNC__ 669d4bb536fSBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" 670e1311b90SBarry Smith int MatDestroy_MPIAIJ(Mat mat) 6711eb62cbbSBarry Smith { 67244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 6731eb62cbbSBarry Smith int ierr; 67483e2fdc7SBarry Smith 6753a40ed3dSBarry Smith PetscFunctionBegin; 67670429bc8SBarry Smith if (--mat->refct > 0) PetscFunctionReturn(0); 67770429bc8SBarry Smith 67870429bc8SBarry Smith if (mat->mapping) { 67970429bc8SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 68070429bc8SBarry Smith } 68170429bc8SBarry Smith if (mat->bmapping) { 68270429bc8SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr); 68370429bc8SBarry Smith } 68461b13de0SBarry Smith if (mat->rmap) { 68561b13de0SBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 68661b13de0SBarry Smith } 68761b13de0SBarry Smith if (mat->cmap) { 68861b13de0SBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 68961b13de0SBarry Smith } 6903a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 691e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",aij->M,aij->N); 692a5a9c739SBarry Smith #endif 6938798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash); CHKERRQ(ierr); 6940452661fSBarry Smith PetscFree(aij->rowners); 69578b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 69678b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 697b1fc9764SSatish Balay #if defined (USE_CTABLE) 698b1fc9764SSatish Balay if (aij->colmap) TableDelete(aij->colmap); 699b1fc9764SSatish Balay #else 7000452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 701b1fc9764SSatish Balay #endif 7020452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 7031eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 704a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 7057a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 7060452661fSBarry Smith PetscFree(aij); 707a5a9c739SBarry Smith PLogObjectDestroy(mat); 7080452661fSBarry Smith PetscHeaderDestroy(mat); 7093a40ed3dSBarry Smith PetscFunctionReturn(0); 7101eb62cbbSBarry Smith } 711ee50ffe9SBarry Smith 7125615d1e5SSatish Balay #undef __FUNC__ 713d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" 714bc5ccf88SSatish Balay int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 7151eb62cbbSBarry Smith { 716416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 717416022c9SBarry Smith int ierr; 718416022c9SBarry Smith 7193a40ed3dSBarry Smith PetscFunctionBegin; 72017699dbbSLois Curfman McInnes if (aij->size == 1) { 721416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 722416022c9SBarry Smith } 723a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 7243a40ed3dSBarry Smith PetscFunctionReturn(0); 725416022c9SBarry Smith } 726416022c9SBarry Smith 7275615d1e5SSatish Balay #undef __FUNC__ 7287b2a1423SBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworSocket" 729bc5ccf88SSatish Balay int MatView_MPIAIJ_ASCIIorDraworSocket(Mat mat,Viewer viewer) 730416022c9SBarry Smith { 73144a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 732dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 73377ed5343SBarry Smith int ierr, format,shift = C->indexshift,rank = aij->rank ; 73477ed5343SBarry Smith int size = aij->size; 735d636dbe3SBarry Smith FILE *fd; 73619bcc07fSBarry Smith ViewerType vtype; 737416022c9SBarry Smith 7383a40ed3dSBarry Smith PetscFunctionBegin; 73919bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 7403f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 741d8467735SBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 7420a198c4cSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 7434e220ebcSLois Curfman McInnes MatInfo info; 7444e220ebcSLois Curfman McInnes int flg; 745a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 74690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 7474e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 74895e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 74977c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 75095e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 7514e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 75295e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 7534e220ebcSLois Curfman McInnes rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 7544e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->A,MAT_LOCAL,&info); 7554e220ebcSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used); 7564e220ebcSLois Curfman McInnes ierr = MatGetInfo(aij->B,MAT_LOCAL,&info); 7574e220ebcSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used); 758a56f8943SBarry Smith fflush(fd); 75977c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 760a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 7613a40ed3dSBarry Smith PetscFunctionReturn(0); 7623a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 7633a40ed3dSBarry Smith PetscFunctionReturn(0); 76408480c60SBarry Smith } 7653f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 76619bcc07fSBarry Smith Draw draw; 76719bcc07fSBarry Smith PetscTruth isnull; 76877ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw); CHKERRQ(ierr); 7693a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 77019bcc07fSBarry Smith } 77119bcc07fSBarry Smith 77217699dbbSLois Curfman McInnes if (size == 1) { 77378b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 7743a40ed3dSBarry Smith } else { 77595373324SBarry Smith /* assemble the entire matrix onto first processor. */ 77695373324SBarry Smith Mat A; 777ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 7782eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 77995373324SBarry Smith Scalar *a; 7802ee70a88SLois Curfman McInnes 78117699dbbSLois Curfman McInnes if (!rank) { 78255843e3eSBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr); 7833a40ed3dSBarry Smith } else { 78455843e3eSBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr); 78595373324SBarry Smith } 786464493b3SBarry Smith PLogObjectParent(mat,A); 787416022c9SBarry Smith 78895373324SBarry Smith /* copy over the A part */ 789ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 7902ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 79195373324SBarry Smith row = aij->rstart; 792dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 79395373324SBarry Smith for ( i=0; i<m; i++ ) { 794416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 79595373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 79695373324SBarry Smith } 7972ee70a88SLois Curfman McInnes aj = Aloc->j; 798dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 79995373324SBarry Smith 80095373324SBarry Smith /* copy over the B part */ 801ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 8022ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 80395373324SBarry Smith row = aij->rstart; 8040452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 805dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 80695373324SBarry Smith for ( i=0; i<m; i++ ) { 807416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 80895373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 80995373324SBarry Smith } 8100452661fSBarry Smith PetscFree(ct); 8116d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8126d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 81355843e3eSBarry Smith /* 81455843e3eSBarry Smith Everyone has to call to draw the matrix since the graphics waits are 81555843e3eSBarry Smith synchronized across all processors that share the Draw object 81655843e3eSBarry Smith */ 8173f1db9ecSBarry Smith if (!rank || PetscTypeCompare(vtype,DRAW_VIEWER)) { 81878b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 81995373324SBarry Smith } 82078b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 82195373324SBarry Smith } 8223a40ed3dSBarry Smith PetscFunctionReturn(0); 8231eb62cbbSBarry Smith } 8241eb62cbbSBarry Smith 8255615d1e5SSatish Balay #undef __FUNC__ 826d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ" 827e1311b90SBarry Smith int MatView_MPIAIJ(Mat mat,Viewer viewer) 828416022c9SBarry Smith { 829416022c9SBarry Smith int ierr; 83019bcc07fSBarry Smith ViewerType vtype; 831416022c9SBarry Smith 8323a40ed3dSBarry Smith PetscFunctionBegin; 83319bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 8343f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER) || PetscTypeCompare(vtype,DRAW_VIEWER) || 8357b2a1423SBarry Smith PetscTypeCompare(vtype,SOCKET_VIEWER)) { 8367b2a1423SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworSocket(mat,viewer); CHKERRQ(ierr); 8373f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 8383a40ed3dSBarry Smith ierr = MatView_MPIAIJ_Binary(mat,viewer);CHKERRQ(ierr); 8395cd90555SBarry Smith } else { 8405cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 841416022c9SBarry Smith } 8423a40ed3dSBarry Smith PetscFunctionReturn(0); 843416022c9SBarry Smith } 844416022c9SBarry Smith 8451eb62cbbSBarry Smith /* 8461eb62cbbSBarry Smith This has to provide several versions. 8471eb62cbbSBarry Smith 8481eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 8491eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 850d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 8511eb62cbbSBarry Smith */ 8525615d1e5SSatish Balay #undef __FUNC__ 8535615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ" 8548f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 855dbb450caSBarry Smith double fshift,int its,Vec xx) 8568a729477SBarry Smith { 85744a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 858d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 859ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 860c16cb8f2SBarry Smith Scalar *b,*x,*xs,*ls,d,*v,sum; 8616abc6512SBarry Smith int ierr,*idx, *diag; 862416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 8638a729477SBarry Smith 8643a40ed3dSBarry Smith PetscFunctionBegin; 8652e8a6d31SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 8662e8a6d31SBarry Smith ierr = VecGetArray(bb,&b); CHKERRQ(ierr); 8672e8a6d31SBarry Smith ierr = VecGetArray(mat->lvec,&ls);CHKERRQ(ierr); 868dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 869dbb450caSBarry Smith ls = ls + shift; 87083e2fdc7SBarry Smith if (!A->diag) {ierr = MatMarkDiag_SeqAIJ(AA); CHKERRQ(ierr);} 871d6dfbf8fSBarry Smith diag = A->diag; 872c16cb8f2SBarry Smith if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 873da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 874f830108cSBarry Smith ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 8752e8a6d31SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 8762e8a6d31SBarry Smith ierr = VecRestoreArray(bb,&b); CHKERRQ(ierr); 8772e8a6d31SBarry Smith ierr = VecRestoreArray(mat->lvec,&ls);CHKERRQ(ierr); 8783a40ed3dSBarry Smith PetscFunctionReturn(0); 879da3a660dSBarry Smith } 8803a40ed3dSBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 8813a40ed3dSBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 882d6dfbf8fSBarry Smith while (its--) { 883d6dfbf8fSBarry Smith /* go down through the rows */ 884d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 885d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 8864d197716SBarry Smith PLogFlops(4*n+3); 887dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 888dbb450caSBarry Smith v = A->a + A->i[i] + shift; 889d6dfbf8fSBarry Smith sum = b[i]; 890d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 891dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 892d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 893dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 894dbb450caSBarry Smith v = B->a + B->i[i] + shift; 895d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 89655a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 897d6dfbf8fSBarry Smith } 898d6dfbf8fSBarry Smith /* come up through the rows */ 899d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 900d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 9014d197716SBarry Smith PLogFlops(4*n+3) 902dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 903dbb450caSBarry Smith v = A->a + A->i[i] + shift; 904d6dfbf8fSBarry Smith sum = b[i]; 905d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 906dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 907d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 908dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 909dbb450caSBarry Smith v = B->a + B->i[i] + shift; 910d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 91155a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 912d6dfbf8fSBarry Smith } 913d6dfbf8fSBarry Smith } 9143a40ed3dSBarry Smith } else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 915da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 916f830108cSBarry Smith ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9172e8a6d31SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 9182e8a6d31SBarry Smith ierr = VecRestoreArray(bb,&b); CHKERRQ(ierr); 9192e8a6d31SBarry Smith ierr = VecRestoreArray(mat->lvec,&ls);CHKERRQ(ierr); 9203a40ed3dSBarry Smith PetscFunctionReturn(0); 921da3a660dSBarry Smith } 9223a40ed3dSBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 9233a40ed3dSBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 924d6dfbf8fSBarry Smith while (its--) { 925d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 926d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 9274d197716SBarry Smith PLogFlops(4*n+3); 928dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 929dbb450caSBarry Smith v = A->a + A->i[i] + shift; 930d6dfbf8fSBarry Smith sum = b[i]; 931d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 932dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 933d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 934dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 935dbb450caSBarry Smith v = B->a + B->i[i] + shift; 936d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 93755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 938d6dfbf8fSBarry Smith } 939d6dfbf8fSBarry Smith } 9403a40ed3dSBarry Smith } else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 941da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 942f830108cSBarry Smith ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr); 9432e8a6d31SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 9442e8a6d31SBarry Smith ierr = VecRestoreArray(bb,&b); CHKERRQ(ierr); 9452e8a6d31SBarry Smith ierr = VecRestoreArray(mat->lvec,&ls);CHKERRQ(ierr); 9463a40ed3dSBarry Smith PetscFunctionReturn(0); 947da3a660dSBarry Smith } 9482e8a6d31SBarry Smith ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 9492e8a6d31SBarry Smith ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr); 950d6dfbf8fSBarry Smith while (its--) { 951d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 952d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 9534d197716SBarry Smith PLogFlops(4*n+3); 954dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 955dbb450caSBarry Smith v = A->a + A->i[i] + shift; 956d6dfbf8fSBarry Smith sum = b[i]; 957d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 958dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 959d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 960dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 961dbb450caSBarry Smith v = B->a + B->i[i] + shift; 962d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 96355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d; 964d6dfbf8fSBarry Smith } 965d6dfbf8fSBarry Smith } 9663a40ed3dSBarry Smith } else { 967a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Parallel SOR not supported"); 968c16cb8f2SBarry Smith } 9692e8a6d31SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 9702e8a6d31SBarry Smith ierr = VecRestoreArray(bb,&b); CHKERRQ(ierr); 9712e8a6d31SBarry Smith ierr = VecRestoreArray(mat->lvec,&ls);CHKERRQ(ierr); 9723a40ed3dSBarry Smith PetscFunctionReturn(0); 9738a729477SBarry Smith } 974a66be287SLois Curfman McInnes 9755615d1e5SSatish Balay #undef __FUNC__ 976d4bb536fSBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" 9778f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info) 978a66be287SLois Curfman McInnes { 979a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 980a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 9814e220ebcSLois Curfman McInnes int ierr; 9824e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 983a66be287SLois Curfman McInnes 9843a40ed3dSBarry Smith PetscFunctionBegin; 9854e220ebcSLois Curfman McInnes info->block_size = 1.0; 9864e220ebcSLois Curfman McInnes ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr); 9874e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 9884e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 9894e220ebcSLois Curfman McInnes ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr); 9904e220ebcSLois Curfman McInnes isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded; 9914e220ebcSLois Curfman McInnes isend[3] += info->memory; isend[4] += info->mallocs; 992a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 9934e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 9944e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 9954e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 9964e220ebcSLois Curfman McInnes info->memory = isend[3]; 9974e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 998a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 999ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm);CHKERRQ(ierr); 10004e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10014e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10024e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10034e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10044e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1005a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 1006ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm);CHKERRQ(ierr); 10074e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 10084e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 10094e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 10104e220ebcSLois Curfman McInnes info->memory = irecv[3]; 10114e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 1012a66be287SLois Curfman McInnes } 10134e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 10144e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10154e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10168d700155SBarry Smith info->rows_global = (double)mat->M; 10178d700155SBarry Smith info->columns_global = (double)mat->N; 10188d700155SBarry Smith info->rows_local = (double)mat->m; 10198d700155SBarry Smith info->columns_local = (double)mat->N; 10204e220ebcSLois Curfman McInnes 10213a40ed3dSBarry Smith PetscFunctionReturn(0); 1022a66be287SLois Curfman McInnes } 1023a66be287SLois Curfman McInnes 10245615d1e5SSatish Balay #undef __FUNC__ 1025d4bb536fSBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" 10268f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op) 1027c74985f6SBarry Smith { 1028c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 10291dee9f54SBarry Smith int ierr; 1030c74985f6SBarry Smith 10313a40ed3dSBarry Smith PetscFunctionBegin; 10326d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 10336d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 10346da5968aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 1035c2653b3dSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 10364787f768SSatish Balay op == MAT_NEW_NONZERO_ALLOCATION_ERR || 10374787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR) { 10381dee9f54SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 10391dee9f54SBarry Smith ierr = MatSetOption(a->B,op);CHKERRQ(ierr); 1040b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 1041aeafbbfcSLois Curfman McInnes a->roworiented = 1; 10421dee9f54SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 10431dee9f54SBarry Smith ierr = MatSetOption(a->B,op);CHKERRQ(ierr); 1044b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 10456da5968aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 10466d4a8577SBarry Smith op == MAT_SYMMETRIC || 10476d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 10486d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 1049981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIAIJ:Option ignored\n"); 10506d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) { 10514b0e389bSBarry Smith a->roworiented = 0; 10521dee9f54SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 10531dee9f54SBarry Smith ierr = MatSetOption(a->B,op);CHKERRQ(ierr); 10542b362799SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 105590f02eecSBarry Smith a->donotstash = 1; 10563a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS){ 10573a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 10583a40ed3dSBarry Smith } else { 10593a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 10603a40ed3dSBarry Smith } 10613a40ed3dSBarry Smith PetscFunctionReturn(0); 1062c74985f6SBarry Smith } 1063c74985f6SBarry Smith 10645615d1e5SSatish Balay #undef __FUNC__ 1065d4bb536fSBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" 10668f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1067c74985f6SBarry Smith { 106844a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 10693a40ed3dSBarry Smith 10703a40ed3dSBarry Smith PetscFunctionBegin; 10710752156aSBarry Smith if (m) *m = mat->M; 10720752156aSBarry Smith if (n) *n = mat->N; 10733a40ed3dSBarry Smith PetscFunctionReturn(0); 1074c74985f6SBarry Smith } 1075c74985f6SBarry Smith 10765615d1e5SSatish Balay #undef __FUNC__ 1077d4bb536fSBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" 10788f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1079c74985f6SBarry Smith { 108044a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 10813a40ed3dSBarry Smith 10823a40ed3dSBarry Smith PetscFunctionBegin; 10830752156aSBarry Smith if (m) *m = mat->m; 1084f830108cSBarry Smith if (n) *n = mat->n; 10853a40ed3dSBarry Smith PetscFunctionReturn(0); 1086c74985f6SBarry Smith } 1087c74985f6SBarry Smith 10885615d1e5SSatish Balay #undef __FUNC__ 1089d4bb536fSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" 10908f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1091c74985f6SBarry Smith { 109244a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 10933a40ed3dSBarry Smith 10943a40ed3dSBarry Smith PetscFunctionBegin; 1095c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 10963a40ed3dSBarry Smith PetscFunctionReturn(0); 1097c74985f6SBarry Smith } 1098c74985f6SBarry Smith 10995615d1e5SSatish Balay #undef __FUNC__ 11005615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ" 11016d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 110239e00950SLois Curfman McInnes { 1103154123eaSLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 110470f0671dSBarry Smith Scalar *vworkA, *vworkB, **pvA, **pvB,*v_p; 1105154123eaSLois Curfman McInnes int i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart; 1106154123eaSLois Curfman McInnes int nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend; 110770f0671dSBarry Smith int *cmap, *idx_p; 110839e00950SLois Curfman McInnes 11093a40ed3dSBarry Smith PetscFunctionBegin; 1110a8c6a408SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Already active"); 11117a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11127a0afa10SBarry Smith 111370f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11147a0afa10SBarry Smith /* 11157a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11167a0afa10SBarry Smith */ 11177a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 1118c16cb8f2SBarry Smith int max = 1,m = mat->m,tmp; 1119c16cb8f2SBarry Smith for ( i=0; i<m; i++ ) { 11207a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11217a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11227a0afa10SBarry Smith } 11233f97c4b0SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc(max*(sizeof(int)+sizeof(Scalar)));CHKPTRQ(mat->rowvalues); 11247a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11257a0afa10SBarry Smith } 11267a0afa10SBarry Smith 1127a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Only local rows") 1128abc0e9e4SLois Curfman McInnes lrow = row - rstart; 112939e00950SLois Curfman McInnes 1130154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1131154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1132154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 1133f830108cSBarry Smith ierr = (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 1134f830108cSBarry Smith ierr = (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1135154123eaSLois Curfman McInnes nztot = nzA + nzB; 1136154123eaSLois Curfman McInnes 113770f0671dSBarry Smith cmap = mat->garray; 1138154123eaSLois Curfman McInnes if (v || idx) { 1139154123eaSLois Curfman McInnes if (nztot) { 1140154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 114170f0671dSBarry Smith int imark = -1; 1142154123eaSLois Curfman McInnes if (v) { 114370f0671dSBarry Smith *v = v_p = mat->rowvalues; 114439e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 114570f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1146154123eaSLois Curfman McInnes else break; 1147154123eaSLois Curfman McInnes } 1148154123eaSLois Curfman McInnes imark = i; 114970f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 115070f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1151154123eaSLois Curfman McInnes } 1152154123eaSLois Curfman McInnes if (idx) { 115370f0671dSBarry Smith *idx = idx_p = mat->rowindices; 115470f0671dSBarry Smith if (imark > -1) { 115570f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 115670f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 115770f0671dSBarry Smith } 115870f0671dSBarry Smith } else { 1159154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 116070f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1161154123eaSLois Curfman McInnes else break; 1162154123eaSLois Curfman McInnes } 1163154123eaSLois Curfman McInnes imark = i; 116470f0671dSBarry Smith } 116570f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 116670f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 116739e00950SLois Curfman McInnes } 11683f97c4b0SBarry Smith } else { 11691ca473b0SSatish Balay if (idx) *idx = 0; 11701ca473b0SSatish Balay if (v) *v = 0; 11711ca473b0SSatish Balay } 1172154123eaSLois Curfman McInnes } 117339e00950SLois Curfman McInnes *nz = nztot; 1174f830108cSBarry Smith ierr = (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 1175f830108cSBarry Smith ierr = (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 11763a40ed3dSBarry Smith PetscFunctionReturn(0); 117739e00950SLois Curfman McInnes } 117839e00950SLois Curfman McInnes 11795615d1e5SSatish Balay #undef __FUNC__ 1180d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" 11816d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 118239e00950SLois Curfman McInnes { 11837a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 11843a40ed3dSBarry Smith 11853a40ed3dSBarry Smith PetscFunctionBegin; 11867a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 1187a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"MatGetRow not called"); 11887a0afa10SBarry Smith } 11897a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 11903a40ed3dSBarry Smith PetscFunctionReturn(0); 119139e00950SLois Curfman McInnes } 119239e00950SLois Curfman McInnes 11935615d1e5SSatish Balay #undef __FUNC__ 11945615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ" 11958f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1196855ac2c5SLois Curfman McInnes { 1197855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1198ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1199416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1200416022c9SBarry Smith double sum = 0.0; 120104ca555eSLois Curfman McInnes Scalar *v; 120204ca555eSLois Curfman McInnes 12033a40ed3dSBarry Smith PetscFunctionBegin; 120417699dbbSLois Curfman McInnes if (aij->size == 1) { 120514183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 120637fa93a5SLois Curfman McInnes } else { 120704ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 120804ca555eSLois Curfman McInnes v = amat->a; 120904ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 12103a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 1211e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 121204ca555eSLois Curfman McInnes #else 121304ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 121404ca555eSLois Curfman McInnes #endif 121504ca555eSLois Curfman McInnes } 121604ca555eSLois Curfman McInnes v = bmat->a; 121704ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 12183a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 1219e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 122004ca555eSLois Curfman McInnes #else 122104ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 122204ca555eSLois Curfman McInnes #endif 122304ca555eSLois Curfman McInnes } 1224ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr); 122504ca555eSLois Curfman McInnes *norm = sqrt(*norm); 12263a40ed3dSBarry Smith } else if (type == NORM_1) { /* max column norm */ 122704ca555eSLois Curfman McInnes double *tmp, *tmp2; 122804ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 1229758f045eSSatish Balay tmp = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp); 1230758f045eSSatish Balay tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2); 1231cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 123204ca555eSLois Curfman McInnes *norm = 0.0; 123304ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 123404ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1235579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 123604ca555eSLois Curfman McInnes } 123704ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 123804ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1239579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 124004ca555eSLois Curfman McInnes } 1241ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr); 124204ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 124304ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 124404ca555eSLois Curfman McInnes } 12450452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 12463a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 1247515d9167SLois Curfman McInnes double ntemp = 0.0; 124804ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1249dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 125004ca555eSLois Curfman McInnes sum = 0.0; 125104ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1252cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 125304ca555eSLois Curfman McInnes } 1254dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 125504ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1256cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 125704ca555eSLois Curfman McInnes } 1258515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 125904ca555eSLois Curfman McInnes } 1260ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);CHKERRQ(ierr); 1261ca161407SBarry Smith } else { 1262a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 126304ca555eSLois Curfman McInnes } 126437fa93a5SLois Curfman McInnes } 12653a40ed3dSBarry Smith PetscFunctionReturn(0); 1266855ac2c5SLois Curfman McInnes } 1267855ac2c5SLois Curfman McInnes 12685615d1e5SSatish Balay #undef __FUNC__ 12695615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ" 12708f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1271b7c46309SBarry Smith { 1272b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1273dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1274416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1275b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 12763a40ed3dSBarry Smith Mat B; 1277b7c46309SBarry Smith Scalar *array; 1278b7c46309SBarry Smith 12793a40ed3dSBarry Smith PetscFunctionBegin; 1280d4bb536fSBarry Smith if (matout == PETSC_NULL && M != N) { 1281a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place"); 1282d4bb536fSBarry Smith } 1283d4bb536fSBarry Smith 1284d4bb536fSBarry Smith ierr = MatCreateMPIAIJ(A->comm,a->n,a->m,N,M,0,PETSC_NULL,0,PETSC_NULL,&B);CHKERRQ(ierr); 1285b7c46309SBarry Smith 1286b7c46309SBarry Smith /* copy over the A part */ 1287ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1288b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1289b7c46309SBarry Smith row = a->rstart; 1290dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1291b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1292416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1293b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1294b7c46309SBarry Smith } 1295b7c46309SBarry Smith aj = Aloc->j; 12964af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1297b7c46309SBarry Smith 1298b7c46309SBarry Smith /* copy over the B part */ 1299ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1300b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1301b7c46309SBarry Smith row = a->rstart; 13020452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1303dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1304b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1305416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1306b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1307b7c46309SBarry Smith } 13080452661fSBarry Smith PetscFree(ct); 13096d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13106d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13113638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13120de55854SLois Curfman McInnes *matout = B; 13130de55854SLois Curfman McInnes } else { 1314f830108cSBarry Smith PetscOps *Abops; 1315f830108cSBarry Smith struct _MatOps *Aops; 1316f830108cSBarry Smith 13170de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13180452661fSBarry Smith PetscFree(a->rowners); 13190de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13200de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 1321b1fc9764SSatish Balay #if defined (USE_CTABLE) 1322b1fc9764SSatish Balay if (a->colmap) TableDelete(a->colmap); 1323b1fc9764SSatish Balay #else 13240452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 1325b1fc9764SSatish Balay #endif 13260452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13270de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1328a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13290452661fSBarry Smith PetscFree(a); 1330f830108cSBarry Smith 1331f830108cSBarry Smith /* 1332f830108cSBarry Smith This is horrible, horrible code. We need to keep the 1333f830108cSBarry Smith A pointers for the bops and ops but copy everything 1334f830108cSBarry Smith else from C. 1335f830108cSBarry Smith */ 1336f830108cSBarry Smith Abops = A->bops; 1337f830108cSBarry Smith Aops = A->ops; 1338f09e8eb9SSatish Balay PetscMemcpy(A,B,sizeof(struct _p_Mat)); 1339f830108cSBarry Smith A->bops = Abops; 1340f830108cSBarry Smith A->ops = Aops; 13410452661fSBarry Smith PetscHeaderDestroy(B); 13420de55854SLois Curfman McInnes } 13433a40ed3dSBarry Smith PetscFunctionReturn(0); 1344b7c46309SBarry Smith } 1345b7c46309SBarry Smith 13465615d1e5SSatish Balay #undef __FUNC__ 13475615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ" 13484b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr) 1349a008b906SSatish Balay { 13504b967eb1SSatish Balay Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 13514b967eb1SSatish Balay Mat a = aij->A, b = aij->B; 1352a008b906SSatish Balay int ierr,s1,s2,s3; 1353a008b906SSatish Balay 13543a40ed3dSBarry Smith PetscFunctionBegin; 13554b967eb1SSatish Balay ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr); 13564b967eb1SSatish Balay if (rr) { 1357e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&s1);CHKERRQ(ierr); 1358a8c6a408SBarry Smith if (s1!=s3) SETERRQ(PETSC_ERR_ARG_SIZ,0,"right vector non-conforming local size"); 13594b967eb1SSatish Balay /* Overlap communication with computation. */ 136043a90d84SBarry Smith ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1361a008b906SSatish Balay } 13624b967eb1SSatish Balay if (ll) { 1363e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&s1);CHKERRQ(ierr); 1364a8c6a408SBarry Smith if (s1!=s2) SETERRQ(PETSC_ERR_ARG_SIZ,0,"left vector non-conforming local size"); 1365f830108cSBarry Smith ierr = (*b->ops->diagonalscale)(b,ll,0); CHKERRQ(ierr); 13664b967eb1SSatish Balay } 13674b967eb1SSatish Balay /* scale the diagonal block */ 1368f830108cSBarry Smith ierr = (*a->ops->diagonalscale)(a,ll,rr); CHKERRQ(ierr); 13694b967eb1SSatish Balay 13704b967eb1SSatish Balay if (rr) { 13714b967eb1SSatish Balay /* Do a scatter end and then right scale the off-diagonal block */ 137243a90d84SBarry Smith ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr); 1373f830108cSBarry Smith ierr = (*b->ops->diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr); 13744b967eb1SSatish Balay } 13754b967eb1SSatish Balay 13763a40ed3dSBarry Smith PetscFunctionReturn(0); 1377a008b906SSatish Balay } 1378a008b906SSatish Balay 1379a008b906SSatish Balay 13805615d1e5SSatish Balay #undef __FUNC__ 1381d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" 13828f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A) 1383682d7d0cSBarry Smith { 1384682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 13853a40ed3dSBarry Smith int ierr; 1386682d7d0cSBarry Smith 13873a40ed3dSBarry Smith PetscFunctionBegin; 13883a40ed3dSBarry Smith if (!a->rank) { 13893a40ed3dSBarry Smith ierr = MatPrintHelp_SeqAIJ(a->A);CHKERRQ(ierr); 13903a40ed3dSBarry Smith } 13913a40ed3dSBarry Smith PetscFunctionReturn(0); 1392682d7d0cSBarry Smith } 1393682d7d0cSBarry Smith 13945615d1e5SSatish Balay #undef __FUNC__ 1395d4bb536fSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" 13968f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs) 13975a838052SSatish Balay { 13983a40ed3dSBarry Smith PetscFunctionBegin; 13995a838052SSatish Balay *bs = 1; 14003a40ed3dSBarry Smith PetscFunctionReturn(0); 14015a838052SSatish Balay } 14025615d1e5SSatish Balay #undef __FUNC__ 1403d4bb536fSBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" 14048f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A) 1405bb5a7306SBarry Smith { 1406bb5a7306SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1407bb5a7306SBarry Smith int ierr; 14083a40ed3dSBarry Smith 14093a40ed3dSBarry Smith PetscFunctionBegin; 1410bb5a7306SBarry Smith ierr = MatSetUnfactored(a->A); CHKERRQ(ierr); 14113a40ed3dSBarry Smith PetscFunctionReturn(0); 1412bb5a7306SBarry Smith } 1413bb5a7306SBarry Smith 1414d4bb536fSBarry Smith #undef __FUNC__ 1415d4bb536fSBarry Smith #define __FUNC__ "MatEqual_MPIAIJ" 1416d4bb536fSBarry Smith int MatEqual_MPIAIJ(Mat A, Mat B, PetscTruth *flag) 1417d4bb536fSBarry Smith { 1418d4bb536fSBarry Smith Mat_MPIAIJ *matB = (Mat_MPIAIJ *) B->data,*matA = (Mat_MPIAIJ *) A->data; 1419d4bb536fSBarry Smith Mat a, b, c, d; 1420d4bb536fSBarry Smith PetscTruth flg; 1421d4bb536fSBarry Smith int ierr; 1422d4bb536fSBarry Smith 14233a40ed3dSBarry Smith PetscFunctionBegin; 1424a8c6a408SBarry Smith if (B->type != MATMPIAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 1425d4bb536fSBarry Smith a = matA->A; b = matA->B; 1426d4bb536fSBarry Smith c = matB->A; d = matB->B; 1427d4bb536fSBarry Smith 1428d4bb536fSBarry Smith ierr = MatEqual(a, c, &flg); CHKERRQ(ierr); 1429d4bb536fSBarry Smith if (flg == PETSC_TRUE) { 1430d4bb536fSBarry Smith ierr = MatEqual(b, d, &flg); CHKERRQ(ierr); 1431d4bb536fSBarry Smith } 1432ca161407SBarry Smith ierr = MPI_Allreduce(&flg, flag, 1, MPI_INT, MPI_LAND, A->comm);CHKERRQ(ierr); 14333a40ed3dSBarry Smith PetscFunctionReturn(0); 1434d4bb536fSBarry Smith } 1435d4bb536fSBarry Smith 1436cb5b572fSBarry Smith #undef __FUNC__ 1437cb5b572fSBarry Smith #define __FUNC__ "MatCopy_MPIAIJ" 1438cb5b572fSBarry Smith int MatCopy_MPIAIJ(Mat A,Mat B,MatStructure str) 1439cb5b572fSBarry Smith { 1440cb5b572fSBarry Smith int ierr; 1441cb5b572fSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data; 1442cb5b572fSBarry Smith Mat_MPIAIJ *b = (Mat_MPIAIJ *)B->data; 1443cb5b572fSBarry Smith 1444cb5b572fSBarry Smith PetscFunctionBegin; 1445cb5b572fSBarry Smith if (str != SAME_NONZERO_PATTERN || B->type != MATMPIAIJ) { 1446cb5b572fSBarry Smith /* because of the column compression in the off-processor part of the matrix a->B, 1447cb5b572fSBarry Smith the number of columns in a->B and b->B may be different, hence we cannot call 1448cb5b572fSBarry Smith the MatCopy() directly on the two parts. If need be, we can provide a more 1449cb5b572fSBarry Smith efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices 1450cb5b572fSBarry Smith then copying the submatrices */ 1451cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1452cb5b572fSBarry Smith } else { 1453cb5b572fSBarry Smith ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 1454cb5b572fSBarry Smith ierr = MatCopy(a->B,b->B,str);CHKERRQ(ierr); 1455cb5b572fSBarry Smith } 1456cb5b572fSBarry Smith PetscFunctionReturn(0); 1457cb5b572fSBarry Smith } 1458cb5b572fSBarry Smith 14592e8a6d31SBarry Smith extern int MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat *); 14602f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 14610a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 14627b2a1423SBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatReuse,Mat **); 14637b2a1423SBarry Smith extern int MatGetSubMatrix_MPIAIJ (Mat ,IS,IS,int,MatReuse,Mat *); 146400e6dbe6SBarry Smith 14658a729477SBarry Smith /* -------------------------------------------------------------------*/ 1466cda55fadSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAIJ, 1467cda55fadSBarry Smith MatGetRow_MPIAIJ, 1468cda55fadSBarry Smith MatRestoreRow_MPIAIJ, 1469cda55fadSBarry Smith MatMult_MPIAIJ, 1470cda55fadSBarry Smith MatMultAdd_MPIAIJ, 1471cda55fadSBarry Smith MatMultTrans_MPIAIJ, 1472cda55fadSBarry Smith MatMultTransAdd_MPIAIJ, 1473cda55fadSBarry Smith 0, 1474cda55fadSBarry Smith 0, 1475cda55fadSBarry Smith 0, 1476cda55fadSBarry Smith 0, 1477cda55fadSBarry Smith 0, 1478cda55fadSBarry Smith 0, 147944a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1480b7c46309SBarry Smith MatTranspose_MPIAIJ, 1481cda55fadSBarry Smith MatGetInfo_MPIAIJ, 1482cda55fadSBarry Smith MatEqual_MPIAIJ, 1483cda55fadSBarry Smith MatGetDiagonal_MPIAIJ, 1484cda55fadSBarry Smith MatDiagonalScale_MPIAIJ, 1485cda55fadSBarry Smith MatNorm_MPIAIJ, 1486cda55fadSBarry Smith MatAssemblyBegin_MPIAIJ, 1487cda55fadSBarry Smith MatAssemblyEnd_MPIAIJ, 14881eb62cbbSBarry Smith 0, 1489cda55fadSBarry Smith MatSetOption_MPIAIJ, 1490cda55fadSBarry Smith MatZeroEntries_MPIAIJ, 1491cda55fadSBarry Smith MatZeroRows_MPIAIJ, 1492cda55fadSBarry Smith 0, 1493cda55fadSBarry Smith 0, 1494cda55fadSBarry Smith 0, 1495cda55fadSBarry Smith 0, 1496cda55fadSBarry Smith MatGetSize_MPIAIJ, 1497cda55fadSBarry Smith MatGetLocalSize_MPIAIJ, 1498cda55fadSBarry Smith MatGetOwnershipRange_MPIAIJ, 1499cda55fadSBarry Smith 0, 1500cda55fadSBarry Smith 0, 1501cda55fadSBarry Smith 0, 1502cda55fadSBarry Smith 0, 15032e8a6d31SBarry Smith MatDuplicate_MPIAIJ, 1504cda55fadSBarry Smith 0, 1505cda55fadSBarry Smith 0, 1506cda55fadSBarry Smith 0, 1507cda55fadSBarry Smith 0, 1508cda55fadSBarry Smith 0, 1509cda55fadSBarry Smith MatGetSubMatrices_MPIAIJ, 1510cda55fadSBarry Smith MatIncreaseOverlap_MPIAIJ, 1511cda55fadSBarry Smith MatGetValues_MPIAIJ, 1512cb5b572fSBarry Smith MatCopy_MPIAIJ, 1513052efed2SBarry Smith MatPrintHelp_MPIAIJ, 1514cda55fadSBarry Smith MatScale_MPIAIJ, 1515cda55fadSBarry Smith 0, 1516cda55fadSBarry Smith 0, 1517cda55fadSBarry Smith 0, 1518cda55fadSBarry Smith MatGetBlockSize_MPIAIJ, 1519cda55fadSBarry Smith 0, 1520cda55fadSBarry Smith 0, 1521cda55fadSBarry Smith 0, 1522cda55fadSBarry Smith 0, 1523cda55fadSBarry Smith MatFDColoringCreate_MPIAIJ, 1524cda55fadSBarry Smith 0, 1525cda55fadSBarry Smith MatSetUnfactored_MPIAIJ, 1526cda55fadSBarry Smith 0, 1527cda55fadSBarry Smith 0, 1528cda55fadSBarry Smith MatGetSubMatrix_MPIAIJ, 1529cda55fadSBarry Smith 0, 1530cda55fadSBarry Smith 0, 1531cda55fadSBarry Smith MatGetMaps_Petsc}; 153236ce4990SBarry Smith 15332e8a6d31SBarry Smith /* ----------------------------------------------------------------------------------------*/ 15342e8a6d31SBarry Smith 1535fb2e594dSBarry Smith EXTERN_C_BEGIN 15362e8a6d31SBarry Smith #undef __FUNC__ 15372e8a6d31SBarry Smith #define __FUNC__ "MatStoreValues_MPIAIJ" 15382e8a6d31SBarry Smith int MatStoreValues_MPIAIJ(Mat mat) 15392e8a6d31SBarry Smith { 15402e8a6d31SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 15412e8a6d31SBarry Smith int ierr; 15422e8a6d31SBarry Smith 15432e8a6d31SBarry Smith PetscFunctionBegin; 15442e8a6d31SBarry Smith ierr = MatStoreValues(aij->A);CHKERRQ(ierr); 15452e8a6d31SBarry Smith ierr = MatStoreValues(aij->B);CHKERRQ(ierr); 15462e8a6d31SBarry Smith PetscFunctionReturn(0); 15472e8a6d31SBarry Smith } 1548fb2e594dSBarry Smith EXTERN_C_END 15492e8a6d31SBarry Smith 1550fb2e594dSBarry Smith EXTERN_C_BEGIN 15512e8a6d31SBarry Smith #undef __FUNC__ 15522e8a6d31SBarry Smith #define __FUNC__ "MatRetrieveValues_MPIAIJ" 15532e8a6d31SBarry Smith int MatRetrieveValues_MPIAIJ(Mat mat) 15542e8a6d31SBarry Smith { 15552e8a6d31SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data; 15562e8a6d31SBarry Smith int ierr; 15572e8a6d31SBarry Smith 15582e8a6d31SBarry Smith PetscFunctionBegin; 15592e8a6d31SBarry Smith ierr = MatRetrieveValues(aij->A);CHKERRQ(ierr); 15602e8a6d31SBarry Smith ierr = MatRetrieveValues(aij->B);CHKERRQ(ierr); 15612e8a6d31SBarry Smith PetscFunctionReturn(0); 15622e8a6d31SBarry Smith } 1563fb2e594dSBarry Smith EXTERN_C_END 15648a729477SBarry Smith 156527508adbSBarry Smith #include "pc.h" 156627508adbSBarry Smith EXTERN_C_BEGIN 15675ef9f2a5SBarry Smith extern int MatGetDiagonalBlock_MPIAIJ(Mat,PetscTruth *,MatReuse,Mat *); 156827508adbSBarry Smith EXTERN_C_END 156927508adbSBarry Smith 15705615d1e5SSatish Balay #undef __FUNC__ 15715615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ" 15721987afe7SBarry Smith /*@C 1573ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 15743a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 15753a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 15763a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 15773a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 15788a729477SBarry Smith 1579db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1580db81eaa0SLois Curfman McInnes 15818a729477SBarry Smith Input Parameters: 1582db81eaa0SLois Curfman McInnes + comm - MPI communicator 15837d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 158492e8d321SLois Curfman McInnes This value should be the same as the local size used in creating the 158592e8d321SLois Curfman McInnes y vector for the matrix-vector product y = Ax. 15861a3896d6SBarry Smith . n - This value should be the same as the local size used in creating the 15871a3896d6SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 15881a3896d6SBarry Smith calculated if N is given) For square matrices n is almost always m. 158960d380a7SBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 159060d380a7SBarry Smith . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 1591af1d9917SSatish Balay . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 1592af1d9917SSatish Balay (same value is used for all local rows) 1593af1d9917SSatish Balay . d_nnz - array containing the number of nonzeros in the various rows of the 1594af1d9917SSatish Balay DIAGONAL portion of the local submatrix (possibly different for each row) 1595af1d9917SSatish Balay or PETSC_NULL, if d_nz is used to specify the nonzero structure. 1596af1d9917SSatish Balay The size of this array is equal to the number of local rows, i.e 'm'. 1597af1d9917SSatish Balay You must leave room for the diagonal entry even if it is zero. 1598af1d9917SSatish Balay . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 1599af1d9917SSatish Balay submatrix (same value is used for all local rows). 1600af1d9917SSatish Balay - o_nnz - array containing the number of nonzeros in the various rows of the 1601af1d9917SSatish Balay OFF-DIAGONAL portion of the local submatrix (possibly different for 1602af1d9917SSatish Balay each row) or PETSC_NULL, if o_nz is used to specify the nonzero 1603af1d9917SSatish Balay structure. The size of this array is equal to the number 1604af1d9917SSatish Balay of local rows, i.e 'm'. 16058a729477SBarry Smith 1606ff756334SLois Curfman McInnes Output Parameter: 160744cd7ae7SLois Curfman McInnes . A - the matrix 16088a729477SBarry Smith 1609b259b22eSLois Curfman McInnes Notes: 1610af1d9917SSatish Balay m,n,M,N parameters specify the size of the matrix, and its partitioning across 1611af1d9917SSatish Balay processors, while d_nz,d_nnz,o_nz,o_nnz parameters specify the approximate 1612af1d9917SSatish Balay storage requirements for this matrix. 1613af1d9917SSatish Balay 1614af1d9917SSatish Balay If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one 1615af1d9917SSatish Balay processor than it must be used on all processors that share the object for 1616af1d9917SSatish Balay that argument. 1617be79a94dSBarry Smith 1618ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1619ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 16200002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 16210002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1622ff756334SLois Curfman McInnes 1623ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1624ff756334SLois Curfman McInnes (possibly both). 1625ff756334SLois Curfman McInnes 1626af1d9917SSatish Balay The parallel matrix is partitioned such that the first m0 rows belong to 1627af1d9917SSatish Balay process 0, the next m1 rows belong to process 1, the next m2 rows belong 1628af1d9917SSatish Balay to process 2 etc.. where m0,m1,m2... are the input parameter 'm'. 1629af1d9917SSatish Balay 1630af1d9917SSatish Balay The DIAGONAL portion of the local submatrix of a processor can be defined 1631af1d9917SSatish Balay as the submatrix which is obtained by extraction the part corresponding 1632af1d9917SSatish Balay to the rows r1-r2 and columns r1-r2 of the global matrix, where r1 is the 1633af1d9917SSatish Balay first row that belongs to the processor, and r2 is the last row belonging 1634af1d9917SSatish Balay to the this processor. This is a square mxm matrix. The remaining portion 1635af1d9917SSatish Balay of the local submatrix (mxN) constitute the OFF-DIAGONAL portion. 1636af1d9917SSatish Balay 1637af1d9917SSatish Balay If o_nnz, d_nnz are specified, then o_nz, and d_nz are ignored. 1638af1d9917SSatish Balay 16395511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 16405511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 16415511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 16425511cfe3SLois Curfman McInnes 16435511cfe3SLois Curfman McInnes Options Database Keys: 1644db81eaa0SLois Curfman McInnes + -mat_aij_no_inode - Do not use inodes 1645db81eaa0SLois Curfman McInnes . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5) 1646db81eaa0SLois Curfman McInnes - -mat_aij_oneindex - Internally use indexing starting at 1 1647db81eaa0SLois Curfman McInnes rather than 0. Note that when calling MatSetValues(), 1648db81eaa0SLois Curfman McInnes the user still MUST index entries starting at 0! 1649494eafd4SBarry Smith . -mat_mpi - use the parallel matrix data structures even on one processor 1650494eafd4SBarry Smith (defaults to using SeqBAIJ format on one processor) 1651494eafd4SBarry Smith . -mat_mpi - use the parallel matrix data structures even on one processor 1652494eafd4SBarry Smith (defaults to using SeqAIJ format on one processor) 1653494eafd4SBarry Smith 16545511cfe3SLois Curfman McInnes 1655af1d9917SSatish Balay Example usage: 1656e0245417SLois Curfman McInnes 1657af1d9917SSatish Balay Consider the following 8x8 matrix with 34 non-zero values, that is 1658af1d9917SSatish Balay assembled across 3 processors. Lets assume that proc0 owns 3 rows, 1659af1d9917SSatish Balay proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown 1660af1d9917SSatish Balay as follows: 16612191d07cSBarry Smith 1662db81eaa0SLois Curfman McInnes .vb 1663af1d9917SSatish Balay 1 2 0 | 0 3 0 | 0 4 1664af1d9917SSatish Balay Proc0 0 5 6 | 7 0 0 | 8 0 1665af1d9917SSatish Balay 9 0 10 | 11 0 0 | 12 0 1666af1d9917SSatish Balay ------------------------------------- 1667af1d9917SSatish Balay 13 0 14 | 15 16 17 | 0 0 1668af1d9917SSatish Balay Proc1 0 18 0 | 19 20 21 | 0 0 1669af1d9917SSatish Balay 0 0 0 | 22 23 0 | 24 0 1670af1d9917SSatish Balay ------------------------------------- 1671af1d9917SSatish Balay Proc2 25 26 27 | 0 0 28 | 29 0 1672af1d9917SSatish Balay 30 0 0 | 31 32 33 | 0 34 1673db81eaa0SLois Curfman McInnes .ve 1674b810aeb4SBarry Smith 1675af1d9917SSatish Balay This can be represented as a collection of submatrices as: 16765511cfe3SLois Curfman McInnes 1677af1d9917SSatish Balay .vb 1678af1d9917SSatish Balay A B C 1679af1d9917SSatish Balay D E F 1680af1d9917SSatish Balay G H I 1681af1d9917SSatish Balay .ve 1682af1d9917SSatish Balay 1683af1d9917SSatish Balay Where the submatrices A,B,C are owned by proc0, D,E,F are 1684af1d9917SSatish Balay owned by proc1, G,H,I are owned by proc2. 1685af1d9917SSatish Balay 1686af1d9917SSatish Balay The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively. 1687af1d9917SSatish Balay The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively. 1688af1d9917SSatish Balay The 'M','N' parameters are 8,8, and have the same values on all procs. 1689af1d9917SSatish Balay 1690af1d9917SSatish Balay The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are 1691af1d9917SSatish Balay submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices 1692af1d9917SSatish Balay corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively. 1693af1d9917SSatish Balay Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL 1694af1d9917SSatish Balay part as SeqAIJ matrices. for eg: proc1 will store [E] as a SeqAIJ 1695af1d9917SSatish Balay matrix, ans [DF] as another SeqAIJ matrix. 1696af1d9917SSatish Balay 1697af1d9917SSatish Balay When d_nz, o_nz parameters are specified, d_nz storage elements are 1698af1d9917SSatish Balay allocated for every row of the local diagonal submatrix, and o_nz 1699af1d9917SSatish Balay storage locations are allocated for every row of the OFF-DIAGONAL submat. 1700af1d9917SSatish Balay One way to choose d_nz and o_nz is to use the max nonzerors per local 1701af1d9917SSatish Balay rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices. 1702af1d9917SSatish Balay In this case, the values of d_nz,o_nz are: 1703af1d9917SSatish Balay .vb 1704af1d9917SSatish Balay proc0 : dnz = 2, o_nz = 2 1705af1d9917SSatish Balay proc1 : dnz = 3, o_nz = 2 1706af1d9917SSatish Balay proc2 : dnz = 1, o_nz = 4 1707af1d9917SSatish Balay .ve 1708af1d9917SSatish Balay We are allocating m*(d_nz+o_nz) storage locations for every proc. This 1709af1d9917SSatish Balay translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10 1710af1d9917SSatish Balay for proc3. i.e we are using 12+15+10=37 storage locations to store 1711af1d9917SSatish Balay 34 values. 1712af1d9917SSatish Balay 1713af1d9917SSatish Balay When d_nnz, o_nnz parameters are specified, the storage is specified 1714af1d9917SSatish Balay for every row, coresponding to both DIAGONAL and OFF-DIAGONAL submatrices. 1715af1d9917SSatish Balay In the above case the values for d_nnz,o_nnz are: 1716af1d9917SSatish Balay .vb 1717af1d9917SSatish Balay proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2] 1718af1d9917SSatish Balay proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1] 1719af1d9917SSatish Balay proc2: d_nnz = [1,1] and o_nnz = [4,4] 1720af1d9917SSatish Balay .ve 1721af1d9917SSatish Balay Here the space allocated is sum of all the above values i.e 34, and 1722af1d9917SSatish Balay hence pre-allocation is perfect. 17233a511b96SLois Curfman McInnes 1724027ccd11SLois Curfman McInnes Level: intermediate 1725027ccd11SLois Curfman McInnes 1726dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1727ff756334SLois Curfman McInnes 1728fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 17298a729477SBarry Smith @*/ 1730e1311b90SBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A) 17318a729477SBarry Smith { 173244cd7ae7SLois Curfman McInnes Mat B; 173344cd7ae7SLois Curfman McInnes Mat_MPIAIJ *b; 1734eac7125bSBarry Smith int ierr, i,size,flag1 = 0, flag2 = 0; 1735416022c9SBarry Smith 17363a40ed3dSBarry Smith PetscFunctionBegin; 17373914022bSBarry Smith MPI_Comm_size(comm,&size); 17387be0e774SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_mpiaij",&flag1); CHKERRQ(ierr); 17397be0e774SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_mpi",&flag2); CHKERRQ(ierr); 17407be0e774SBarry Smith if (!flag1 && !flag2 && size == 1) { 17413914022bSBarry Smith if (M == PETSC_DECIDE) M = m; 17423914022bSBarry Smith if (N == PETSC_DECIDE) N = n; 17433914022bSBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,d_nz,d_nnz,A); CHKERRQ(ierr); 17443a40ed3dSBarry Smith PetscFunctionReturn(0); 17453914022bSBarry Smith } 17463914022bSBarry Smith 174744cd7ae7SLois Curfman McInnes *A = 0; 17483f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIAIJ,"Mat",comm,MatDestroy,MatView); 174944cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 175044cd7ae7SLois Curfman McInnes B->data = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b); 175144cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_MPIAIJ)); 1752cda55fadSBarry Smith PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)); 1753e1311b90SBarry Smith B->ops->destroy = MatDestroy_MPIAIJ; 1754e1311b90SBarry Smith B->ops->view = MatView_MPIAIJ; 175544cd7ae7SLois Curfman McInnes B->factor = 0; 175644cd7ae7SLois Curfman McInnes B->assembled = PETSC_FALSE; 175790f02eecSBarry Smith B->mapping = 0; 1758d6dfbf8fSBarry Smith 175947794344SBarry Smith B->insertmode = NOT_SET_VALUES; 17609eb4d147SSatish Balay b->size = size; 176144cd7ae7SLois Curfman McInnes MPI_Comm_rank(comm,&b->rank); 17621eb62cbbSBarry Smith 17633a40ed3dSBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) { 1764a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz"); 17653a40ed3dSBarry Smith } 17661987afe7SBarry Smith 1767eac7125bSBarry Smith ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr); 1768eac7125bSBarry Smith ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr); 176944cd7ae7SLois Curfman McInnes b->m = m; B->m = m; 177044cd7ae7SLois Curfman McInnes b->n = n; B->n = n; 177144cd7ae7SLois Curfman McInnes b->N = N; B->N = N; 177244cd7ae7SLois Curfman McInnes b->M = M; B->M = M; 17731eb62cbbSBarry Smith 1774c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 1775c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 1776253a16ddSBarry Smith ierr = MapCreateMPI(comm,m,M,&B->rmap);CHKERRQ(ierr); 1777253a16ddSBarry Smith ierr = MapCreateMPI(comm,n,N,&B->cmap);CHKERRQ(ierr); 1778c7fcc2eaSBarry Smith 17791eb62cbbSBarry Smith /* build local table of row and column ownerships */ 178044cd7ae7SLois Curfman McInnes b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners); 1781f09e8eb9SSatish Balay PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1782603f58a4SSatish Balay b->cowners = b->rowners + b->size + 2; 1783ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 178444cd7ae7SLois Curfman McInnes b->rowners[0] = 0; 178544cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 178644cd7ae7SLois Curfman McInnes b->rowners[i] += b->rowners[i-1]; 17878a729477SBarry Smith } 178844cd7ae7SLois Curfman McInnes b->rstart = b->rowners[b->rank]; 178944cd7ae7SLois Curfman McInnes b->rend = b->rowners[b->rank+1]; 1790ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 179144cd7ae7SLois Curfman McInnes b->cowners[0] = 0; 179244cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 179344cd7ae7SLois Curfman McInnes b->cowners[i] += b->cowners[i-1]; 17948a729477SBarry Smith } 179544cd7ae7SLois Curfman McInnes b->cstart = b->cowners[b->rank]; 179644cd7ae7SLois Curfman McInnes b->cend = b->cowners[b->rank+1]; 17978a729477SBarry Smith 17985ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 1799029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr); 180044cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->A); 18017b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 1802029af93fSBarry Smith ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr); 180344cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->B); 18048a729477SBarry Smith 18051eb62cbbSBarry Smith /* build cache for off array entries formed */ 18068798bf22SSatish Balay ierr = MatStashCreate_Private(comm,1,&B->stash); CHKERRQ(ierr); 180790f02eecSBarry Smith b->donotstash = 0; 180844cd7ae7SLois Curfman McInnes b->colmap = 0; 180944cd7ae7SLois Curfman McInnes b->garray = 0; 181044cd7ae7SLois Curfman McInnes b->roworiented = 1; 18118a729477SBarry Smith 18121eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 181344cd7ae7SLois Curfman McInnes b->lvec = 0; 181444cd7ae7SLois Curfman McInnes b->Mvctx = 0; 18158a729477SBarry Smith 18167a0afa10SBarry Smith /* stuff for MatGetRow() */ 181744cd7ae7SLois Curfman McInnes b->rowindices = 0; 181844cd7ae7SLois Curfman McInnes b->rowvalues = 0; 181944cd7ae7SLois Curfman McInnes b->getrowactive = PETSC_FALSE; 18207a0afa10SBarry Smith 18212e8a6d31SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 18222e8a6d31SBarry Smith "MatStoreValues_MPIAIJ", 18232e8a6d31SBarry Smith (void*)MatStoreValues_MPIAIJ);CHKERRQ(ierr); 18242e8a6d31SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 18252e8a6d31SBarry Smith "MatRetrieveValues_MPIAIJ", 18262e8a6d31SBarry Smith (void*)MatRetrieveValues_MPIAIJ);CHKERRQ(ierr); 18275ef9f2a5SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetDiagonalBlock_C", 18285ef9f2a5SBarry Smith "MatGetDiagonalBlock_MPIAIJ", 18295ef9f2a5SBarry Smith (void*)MatGetDiagonalBlock_MPIAIJ);CHKERRQ(ierr); 183044cd7ae7SLois Curfman McInnes *A = B; 18313a40ed3dSBarry Smith PetscFunctionReturn(0); 1832d6dfbf8fSBarry Smith } 1833c74985f6SBarry Smith 18345615d1e5SSatish Balay #undef __FUNC__ 18352e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_MPIAIJ" 18362e8a6d31SBarry Smith int MatDuplicate_MPIAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat) 1837d6dfbf8fSBarry Smith { 1838d6dfbf8fSBarry Smith Mat mat; 1839416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 1840a1b97e82SLois Curfman McInnes int ierr, len=0, flg; 1841d6dfbf8fSBarry Smith 18423a40ed3dSBarry Smith PetscFunctionBegin; 1843416022c9SBarry Smith *newmat = 0; 18443f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIAIJ,"Mat",matin->comm,MatDestroy,MatView); 1845a5a9c739SBarry Smith PLogObjectCreate(mat); 18460452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1847cda55fadSBarry Smith PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps)); 1848e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIAIJ; 1849e1311b90SBarry Smith mat->ops->view = MatView_MPIAIJ; 1850d6dfbf8fSBarry Smith mat->factor = matin->factor; 1851c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1852e7641de0SSatish Balay mat->insertmode = NOT_SET_VALUES; 1853d6dfbf8fSBarry Smith 185444cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 185544cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 185644cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 185744cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 1858d6dfbf8fSBarry Smith 1859416022c9SBarry Smith a->rstart = oldmat->rstart; 1860416022c9SBarry Smith a->rend = oldmat->rend; 1861416022c9SBarry Smith a->cstart = oldmat->cstart; 1862416022c9SBarry Smith a->cend = oldmat->cend; 186317699dbbSLois Curfman McInnes a->size = oldmat->size; 186417699dbbSLois Curfman McInnes a->rank = oldmat->rank; 1865e7641de0SSatish Balay a->donotstash = oldmat->donotstash; 1866e7641de0SSatish Balay a->roworiented = oldmat->roworiented; 1867e7641de0SSatish Balay a->rowindices = 0; 1868bcd2baecSBarry Smith a->rowvalues = 0; 1869bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1870d6dfbf8fSBarry Smith 1871603f58a4SSatish Balay a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1872f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ)); 1873603f58a4SSatish Balay a->cowners = a->rowners + a->size + 2; 1874603f58a4SSatish Balay PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int)); 18758798bf22SSatish Balay ierr = MatStashCreate_Private(matin->comm,1,&mat->stash); CHKERRQ(ierr); 18762ee70a88SLois Curfman McInnes if (oldmat->colmap) { 1877b1fc9764SSatish Balay #if defined (USE_CTABLE) 1878fa46199cSSatish Balay ierr = TableCreateCopy(oldmat->colmap,&a->colmap); CHKERRQ(ierr); 1879b1fc9764SSatish Balay #else 18800452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1881416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1882416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1883b1fc9764SSatish Balay #endif 1884416022c9SBarry Smith } else a->colmap = 0; 18853f41c07dSBarry Smith if (oldmat->garray) { 18863f41c07dSBarry Smith len = ((Mat_SeqAIJ *) (oldmat->B->data))->n; 18873f41c07dSBarry Smith a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray); 1888464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 18893f41c07dSBarry Smith if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1890416022c9SBarry Smith } else a->garray = 0; 1891d6dfbf8fSBarry Smith 1892416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1893416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1894a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1895416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 18962e8a6d31SBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A); CHKERRQ(ierr); 1897416022c9SBarry Smith PLogObjectParent(mat,a->A); 18982e8a6d31SBarry Smith ierr = MatDuplicate(oldmat->B,cpvalues,&a->B); CHKERRQ(ierr); 1899416022c9SBarry Smith PLogObjectParent(mat,a->B); 19005dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 190125cdf11fSBarry Smith if (flg) { 1902682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1903682d7d0cSBarry Smith } 190427508adbSBarry Smith ierr = FListDuplicate(matin->qlist,&mat->qlist);CHKERRQ(ierr); 19058a729477SBarry Smith *newmat = mat; 19063a40ed3dSBarry Smith PetscFunctionReturn(0); 19078a729477SBarry Smith } 1908416022c9SBarry Smith 190977c4ece6SBarry Smith #include "sys.h" 1910416022c9SBarry Smith 19115615d1e5SSatish Balay #undef __FUNC__ 19125615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ" 191319bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat) 1914416022c9SBarry Smith { 1915d65a2f8fSBarry Smith Mat A; 1916d65a2f8fSBarry Smith Scalar *vals,*svals; 191719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 1918416022c9SBarry Smith MPI_Status status; 19193a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 192017699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1921d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 1922b362ba68SBarry Smith int tag = ((PetscObject)viewer)->tag,cend,cstart,n; 1923416022c9SBarry Smith 19243a40ed3dSBarry Smith PetscFunctionBegin; 192517699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 192617699dbbSLois Curfman McInnes if (!rank) { 192790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 19280752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr); 1929a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 1930d64ed03dSBarry Smith if (header[3] < 0) { 1931a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix in special format on disk, cannot load as MPIAIJ"); 1932d64ed03dSBarry Smith } 19336c5fab8fSBarry Smith } 19346c5fab8fSBarry Smith 1935ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 1936416022c9SBarry Smith M = header[1]; N = header[2]; 1937416022c9SBarry Smith /* determine ownership of all rows */ 193817699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 19390452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1940ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1941416022c9SBarry Smith rowners[0] = 0; 194217699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1943416022c9SBarry Smith rowners[i] += rowners[i-1]; 1944416022c9SBarry Smith } 194517699dbbSLois Curfman McInnes rstart = rowners[rank]; 194617699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1947416022c9SBarry Smith 1948416022c9SBarry Smith /* distribute row lengths to all processors */ 19490452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1950416022c9SBarry Smith offlens = ourlens + (rend-rstart); 195117699dbbSLois Curfman McInnes if (!rank) { 19520452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 19530752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 19540452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 195517699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1956ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 19570452661fSBarry Smith PetscFree(sndcounts); 19583a40ed3dSBarry Smith } else { 1959ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 1960416022c9SBarry Smith } 1961416022c9SBarry Smith 196217699dbbSLois Curfman McInnes if (!rank) { 1963416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 19640452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1965cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 196617699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1967416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1968416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1969416022c9SBarry Smith } 1970416022c9SBarry Smith } 19710452661fSBarry Smith PetscFree(rowlengths); 1972416022c9SBarry Smith 1973416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1974416022c9SBarry Smith maxnz = 0; 197517699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 19760452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1977416022c9SBarry Smith } 19780452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1979416022c9SBarry Smith 1980416022c9SBarry Smith /* read in my part of the matrix column indices */ 1981416022c9SBarry Smith nz = procsnz[0]; 19820452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 19830752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr); 1984d65a2f8fSBarry Smith 1985d65a2f8fSBarry Smith /* read in every one elses and ship off */ 198617699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1987d65a2f8fSBarry Smith nz = procsnz[i]; 19880752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 1989ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 1990d65a2f8fSBarry Smith } 19910452661fSBarry Smith PetscFree(cols); 19923a40ed3dSBarry Smith } else { 1993416022c9SBarry Smith /* determine buffer space needed for message */ 1994416022c9SBarry Smith nz = 0; 1995416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1996416022c9SBarry Smith nz += ourlens[i]; 1997416022c9SBarry Smith } 19980452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1999416022c9SBarry Smith 2000416022c9SBarry Smith /* receive message of column indices*/ 2001ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 2002ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 2003a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 2004416022c9SBarry Smith } 2005416022c9SBarry Smith 2006b362ba68SBarry Smith /* determine column ownership if matrix is not square */ 2007b362ba68SBarry Smith if (N != M) { 2008b362ba68SBarry Smith n = N/size + ((N % size) > rank); 2009b362ba68SBarry Smith ierr = MPI_Scan(&n,&cend,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2010b362ba68SBarry Smith cstart = cend - n; 2011b362ba68SBarry Smith } else { 2012b362ba68SBarry Smith cstart = rstart; 2013b362ba68SBarry Smith cend = rend; 2014fb2e594dSBarry Smith n = cend - cstart; 2015b362ba68SBarry Smith } 2016b362ba68SBarry Smith 2017416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 2018cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 2019416022c9SBarry Smith jj = 0; 2020416022c9SBarry Smith for ( i=0; i<m; i++ ) { 2021416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 2022b362ba68SBarry Smith if (mycols[jj] < cstart || mycols[jj] >= cend) offlens[i]++; 2023416022c9SBarry Smith jj++; 2024416022c9SBarry Smith } 2025416022c9SBarry Smith } 2026d65a2f8fSBarry Smith 2027d65a2f8fSBarry Smith /* create our matrix */ 2028416022c9SBarry Smith for ( i=0; i<m; i++ ) { 2029416022c9SBarry Smith ourlens[i] -= offlens[i]; 2030416022c9SBarry Smith } 2031b362ba68SBarry Smith ierr = MatCreateMPIAIJ(comm,m,n,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 2032d65a2f8fSBarry Smith A = *newmat; 2033fb2e594dSBarry Smith ierr = MatSetOption(A,MAT_COLUMNS_SORTED); CHKERRQ(ierr); 2034d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 2035d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 2036d65a2f8fSBarry Smith } 2037416022c9SBarry Smith 203817699dbbSLois Curfman McInnes if (!rank) { 20390452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 2040416022c9SBarry Smith 2041416022c9SBarry Smith /* read in my part of the matrix numerical values */ 2042416022c9SBarry Smith nz = procsnz[0]; 20430752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 2044d65a2f8fSBarry Smith 2045d65a2f8fSBarry Smith /* insert into matrix */ 2046d65a2f8fSBarry Smith jj = rstart; 2047d65a2f8fSBarry Smith smycols = mycols; 2048d65a2f8fSBarry Smith svals = vals; 2049d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 2050d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 2051d65a2f8fSBarry Smith smycols += ourlens[i]; 2052d65a2f8fSBarry Smith svals += ourlens[i]; 2053d65a2f8fSBarry Smith jj++; 2054416022c9SBarry Smith } 2055416022c9SBarry Smith 2056d65a2f8fSBarry Smith /* read in other processors and ship out */ 205717699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 2058416022c9SBarry Smith nz = procsnz[i]; 20590752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 2060ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 2061416022c9SBarry Smith } 20620452661fSBarry Smith PetscFree(procsnz); 20633a40ed3dSBarry Smith } else { 2064d65a2f8fSBarry Smith /* receive numeric values */ 20650452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 2066416022c9SBarry Smith 2067d65a2f8fSBarry Smith /* receive message of values*/ 2068ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 2069ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 2070a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 2071d65a2f8fSBarry Smith 2072d65a2f8fSBarry Smith /* insert into matrix */ 2073d65a2f8fSBarry Smith jj = rstart; 2074d65a2f8fSBarry Smith smycols = mycols; 2075d65a2f8fSBarry Smith svals = vals; 2076d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 2077d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 2078d65a2f8fSBarry Smith smycols += ourlens[i]; 2079d65a2f8fSBarry Smith svals += ourlens[i]; 2080d65a2f8fSBarry Smith jj++; 2081d65a2f8fSBarry Smith } 2082d65a2f8fSBarry Smith } 20830452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 2084d65a2f8fSBarry Smith 20856d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 20866d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 20873a40ed3dSBarry Smith PetscFunctionReturn(0); 2088416022c9SBarry Smith } 2089a0ff6018SBarry Smith 209029da9460SBarry Smith #undef __FUNC__ 209129da9460SBarry Smith #define __FUNC__ "MatGetSubMatrix_MPIAIJ" 2092a0ff6018SBarry Smith /* 209329da9460SBarry Smith Not great since it makes two copies of the submatrix, first an SeqAIJ 209429da9460SBarry Smith in local and then by concatenating the local matrices the end result. 209529da9460SBarry Smith Writing it directly would be much like MatGetSubMatrices_MPIAIJ() 2096a0ff6018SBarry Smith */ 20977b2a1423SBarry Smith int MatGetSubMatrix_MPIAIJ(Mat mat,IS isrow,IS iscol,int csize,MatReuse call,Mat *newmat) 2098a0ff6018SBarry Smith { 209900e6dbe6SBarry Smith int ierr, i, m,n,rstart,row,rend,nz,*cwork,size,rank,j; 2100fee21e36SBarry Smith Mat *local,M, Mreuse; 210100e6dbe6SBarry Smith Scalar *vwork,*aa; 210200e6dbe6SBarry Smith MPI_Comm comm = mat->comm; 210300e6dbe6SBarry Smith Mat_SeqAIJ *aij; 210400e6dbe6SBarry Smith int *ii, *jj,nlocal,*dlens,*olens,dlen,olen,jend; 2105a0ff6018SBarry Smith 2106a0ff6018SBarry Smith PetscFunctionBegin; 210700e6dbe6SBarry Smith MPI_Comm_rank(comm,&rank); 210800e6dbe6SBarry Smith MPI_Comm_size(comm,&size); 210900e6dbe6SBarry Smith 2110fee21e36SBarry Smith if (call == MAT_REUSE_MATRIX) { 2111fee21e36SBarry Smith ierr = PetscObjectQuery((PetscObject)*newmat,"SubMatrix",(PetscObject *)&Mreuse);CHKERRQ(ierr); 2112fee21e36SBarry Smith if (!Mreuse) SETERRQ(1,1,"Submatrix passed in was not used before, cannot reuse"); 2113fee21e36SBarry Smith local = &Mreuse; 2114fee21e36SBarry Smith ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_REUSE_MATRIX,&local);CHKERRQ(ierr); 2115fee21e36SBarry Smith } else { 2116a0ff6018SBarry Smith ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr); 2117fee21e36SBarry Smith Mreuse = *local; 2118fee21e36SBarry Smith PetscFree(local); 2119fee21e36SBarry Smith } 2120a0ff6018SBarry Smith 2121a0ff6018SBarry Smith /* 2122a0ff6018SBarry Smith m - number of local rows 2123a0ff6018SBarry Smith n - number of columns (same on all processors) 2124a0ff6018SBarry Smith rstart - first row in new global matrix generated 2125a0ff6018SBarry Smith */ 2126fee21e36SBarry Smith ierr = MatGetSize(Mreuse,&m,&n);CHKERRQ(ierr); 2127a0ff6018SBarry Smith if (call == MAT_INITIAL_MATRIX) { 2128fee21e36SBarry Smith aij = (Mat_SeqAIJ *) (Mreuse)->data; 2129a8c6a408SBarry Smith if (aij->indexshift) SETERRQ(PETSC_ERR_SUP,1,"No support for index shifted matrix"); 213000e6dbe6SBarry Smith ii = aij->i; 213100e6dbe6SBarry Smith jj = aij->j; 213200e6dbe6SBarry Smith 2133a0ff6018SBarry Smith /* 213400e6dbe6SBarry Smith Determine the number of non-zeros in the diagonal and off-diagonal 213500e6dbe6SBarry Smith portions of the matrix in order to do correct preallocation 2136a0ff6018SBarry Smith */ 213700e6dbe6SBarry Smith 213800e6dbe6SBarry Smith /* first get start and end of "diagonal" columns */ 21396a6a5d1dSBarry Smith if (csize == PETSC_DECIDE) { 214000e6dbe6SBarry Smith nlocal = n/size + ((n % size) > rank); 21416a6a5d1dSBarry Smith } else { 21426a6a5d1dSBarry Smith nlocal = csize; 21436a6a5d1dSBarry Smith } 2144ca161407SBarry Smith ierr = MPI_Scan(&nlocal,&rend,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 214500e6dbe6SBarry Smith rstart = rend - nlocal; 21466a6a5d1dSBarry Smith if (rank == size - 1 && rend != n) { 21476a6a5d1dSBarry Smith SETERRQ(1,1,"Local column sizes do not add up to total number of columns"); 21486a6a5d1dSBarry Smith } 214900e6dbe6SBarry Smith 215000e6dbe6SBarry Smith /* next, compute all the lengths */ 215100e6dbe6SBarry Smith dlens = (int *) PetscMalloc( (2*m+1)*sizeof(int) );CHKPTRQ(dlens); 215200e6dbe6SBarry Smith olens = dlens + m; 215300e6dbe6SBarry Smith for ( i=0; i<m; i++ ) { 215400e6dbe6SBarry Smith jend = ii[i+1] - ii[i]; 215500e6dbe6SBarry Smith olen = 0; 215600e6dbe6SBarry Smith dlen = 0; 215700e6dbe6SBarry Smith for ( j=0; j<jend; j++ ) { 215800e6dbe6SBarry Smith if ( *jj < rstart || *jj >= rend) olen++; 215900e6dbe6SBarry Smith else dlen++; 216000e6dbe6SBarry Smith jj++; 216100e6dbe6SBarry Smith } 216200e6dbe6SBarry Smith olens[i] = olen; 216300e6dbe6SBarry Smith dlens[i] = dlen; 216400e6dbe6SBarry Smith } 21652d207970SBarry Smith ierr = MatCreateMPIAIJ(comm,m,nlocal,PETSC_DECIDE,n,0,dlens,0,olens,&M);CHKERRQ(ierr); 216600e6dbe6SBarry Smith PetscFree(dlens); 2167a0ff6018SBarry Smith } else { 2168a0ff6018SBarry Smith int ml,nl; 2169a0ff6018SBarry Smith 2170a0ff6018SBarry Smith M = *newmat; 2171a0ff6018SBarry Smith ierr = MatGetLocalSize(M,&ml,&nl);CHKERRQ(ierr); 2172a8c6a408SBarry Smith if (ml != m) SETERRQ(PETSC_ERR_ARG_SIZ,1,"Previous matrix must be same size/layout as request"); 2173a0ff6018SBarry Smith ierr = MatZeroEntries(M);CHKERRQ(ierr); 2174c48de900SBarry Smith /* 2175c48de900SBarry Smith The next two lines are needed so we may call MatSetValues_MPIAIJ() below directly, 2176c48de900SBarry Smith rather than the slower MatSetValues(). 2177c48de900SBarry Smith */ 2178c48de900SBarry Smith M->was_assembled = PETSC_TRUE; 2179c48de900SBarry Smith M->assembled = PETSC_FALSE; 2180a0ff6018SBarry Smith } 2181a0ff6018SBarry Smith ierr = MatGetOwnershipRange(M,&rstart,&rend); CHKERRQ(ierr); 2182fee21e36SBarry Smith aij = (Mat_SeqAIJ *) (Mreuse)->data; 2183a8c6a408SBarry Smith if (aij->indexshift) SETERRQ(PETSC_ERR_SUP,1,"No support for index shifted matrix"); 218400e6dbe6SBarry Smith ii = aij->i; 218500e6dbe6SBarry Smith jj = aij->j; 218600e6dbe6SBarry Smith aa = aij->a; 2187a0ff6018SBarry Smith for (i=0; i<m; i++) { 2188a0ff6018SBarry Smith row = rstart + i; 218900e6dbe6SBarry Smith nz = ii[i+1] - ii[i]; 219000e6dbe6SBarry Smith cwork = jj; jj += nz; 219100e6dbe6SBarry Smith vwork = aa; aa += nz; 21928c638d02SBarry Smith ierr = MatSetValues_MPIAIJ(M,1,&row,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 2193a0ff6018SBarry Smith } 2194a0ff6018SBarry Smith 2195a0ff6018SBarry Smith ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2196a0ff6018SBarry Smith ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2197a0ff6018SBarry Smith *newmat = M; 2198fee21e36SBarry Smith 2199fee21e36SBarry Smith /* save submatrix used in processor for next request */ 2200fee21e36SBarry Smith if (call == MAT_INITIAL_MATRIX) { 2201fee21e36SBarry Smith ierr = PetscObjectCompose((PetscObject)M,"SubMatrix",(PetscObject)Mreuse);CHKERRQ(ierr); 2202fee21e36SBarry Smith ierr = PetscObjectDereference((PetscObject)Mreuse);CHKERRQ(ierr); 2203fee21e36SBarry Smith } 2204fee21e36SBarry Smith 2205a0ff6018SBarry Smith PetscFunctionReturn(0); 2206a0ff6018SBarry Smith } 2207