1cb512458SBarry Smith #ifndef lint 2*f5eb4b81SSatish Balay static char vcid[] = "$Id: mpiaij.c,v 1.140 1996/04/11 14:07:38 curfman Exp balay $"; 3cb512458SBarry Smith #endif 48a729477SBarry Smith 51eb62cbbSBarry Smith #include "mpiaij.h" 6*f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 7*f5eb4b81SSatish Balay #include "vec/inline/spops.h" 88a729477SBarry Smith 99e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column 109e25ed09SBarry Smith number to the local number in the off-diagonal part of the local 119e25ed09SBarry Smith storage of the matrix. This is done in a non scable way since the 129e25ed09SBarry Smith length of colmap equals the global matrix length. 139e25ed09SBarry Smith */ 14b7c46309SBarry Smith static int CreateColmap_Private(Mat mat) 159e25ed09SBarry Smith { 1644a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 17ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data; 18416022c9SBarry Smith int n = B->n,i,shift = B->indexshift; 19dbb450caSBarry Smith 200452661fSBarry Smith aij->colmap = (int *) PetscMalloc(aij->N*sizeof(int));CHKPTRQ(aij->colmap); 21464493b3SBarry Smith PLogObjectMemory(mat,aij->N*sizeof(int)); 22cddf8d76SBarry Smith PetscMemzero(aij->colmap,aij->N*sizeof(int)); 23dbb450caSBarry Smith for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i-shift; 249e25ed09SBarry Smith return 0; 259e25ed09SBarry Smith } 269e25ed09SBarry Smith 272493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat); 282493cbb0SBarry Smith 2951c98ccdSLois Curfman McInnes static int MatGetReordering_MPIAIJ(Mat mat,MatOrdering type,IS *rperm,IS *cperm) 30299609e3SLois Curfman McInnes { 31299609e3SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 32299609e3SLois Curfman McInnes int ierr; 3317699dbbSLois Curfman McInnes if (aij->size == 1) { 3451c98ccdSLois Curfman McInnes ierr = MatGetReordering(aij->A,type,rperm,cperm); CHKERRQ(ierr); 3548d91487SBarry Smith } else SETERRQ(1,"MatGetReordering_MPIAIJ:not supported in parallel"); 36299609e3SLois Curfman McInnes return 0; 37299609e3SLois Curfman McInnes } 38299609e3SLois Curfman McInnes 394b0e389bSBarry Smith static int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 408a729477SBarry Smith { 4144a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 42dbb450caSBarry Smith Mat_SeqAIJ *C = (Mat_SeqAIJ*) aij->A->data; 434b0e389bSBarry Smith Scalar value; 441eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 451eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 464b0e389bSBarry Smith int shift = C->indexshift,roworiented = aij->roworiented; 478a729477SBarry Smith 4864119d58SLois Curfman McInnes if (aij->insertmode != NOT_SET_VALUES && aij->insertmode != addv) { 4948d91487SBarry Smith SETERRQ(1,"MatSetValues_MPIAIJ:Cannot mix inserts and adds"); 508a729477SBarry Smith } 511eb62cbbSBarry Smith aij->insertmode = addv; 528a729477SBarry Smith for ( i=0; i<m; i++ ) { 534b0e389bSBarry Smith if (im[i] < 0) SETERRQ(1,"MatSetValues_MPIAIJ:Negative row"); 544b0e389bSBarry Smith if (im[i] >= aij->M) SETERRQ(1,"MatSetValues_MPIAIJ:Row too large"); 554b0e389bSBarry Smith if (im[i] >= rstart && im[i] < rend) { 564b0e389bSBarry Smith row = im[i] - rstart; 571eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 584b0e389bSBarry Smith if (in[j] < 0) SETERRQ(1,"MatSetValues_MPIAIJ:Negative column"); 594b0e389bSBarry Smith if (in[j] >= aij->N) SETERRQ(1,"MatSetValues_MPIAIJ:Col too large"); 604b0e389bSBarry Smith if (in[j] >= cstart && in[j] < cend){ 614b0e389bSBarry Smith col = in[j] - cstart; 624b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 634b0e389bSBarry Smith ierr = MatSetValues(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); 641eb62cbbSBarry Smith } 651eb62cbbSBarry Smith else { 66227d817aSBarry Smith if (mat->was_assembled) { 67b7c46309SBarry Smith if (!aij->colmap) {ierr = CreateColmap_Private(mat);CHKERRQ(ierr);} 684b0e389bSBarry Smith col = aij->colmap[in[j]] + shift; 69ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 702493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 714b0e389bSBarry Smith col = in[j]; 72d6dfbf8fSBarry Smith } 739e25ed09SBarry Smith } 744b0e389bSBarry Smith else col = in[j]; 754b0e389bSBarry Smith if (roworiented) value = v[i*n+j]; else value = v[i+j*m]; 764b0e389bSBarry Smith ierr = MatSetValues(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); 771eb62cbbSBarry Smith } 781eb62cbbSBarry Smith } 791eb62cbbSBarry Smith } 801eb62cbbSBarry Smith else { 814b0e389bSBarry Smith if (roworiented) { 824b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr); 834b0e389bSBarry Smith } 844b0e389bSBarry Smith else { 854b0e389bSBarry Smith row = im[i]; 864b0e389bSBarry Smith for ( j=0; j<n; j++ ) { 874b0e389bSBarry Smith ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr); 884b0e389bSBarry Smith } 894b0e389bSBarry Smith } 901eb62cbbSBarry Smith } 918a729477SBarry Smith } 928a729477SBarry Smith return 0; 938a729477SBarry Smith } 948a729477SBarry Smith 95b49de8d1SLois Curfman McInnes static int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 96b49de8d1SLois Curfman McInnes { 97b49de8d1SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 98b49de8d1SLois Curfman McInnes Mat_SeqAIJ *C = (Mat_SeqAIJ*) aij->A->data; 99b49de8d1SLois Curfman McInnes int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 100b49de8d1SLois Curfman McInnes int cstart = aij->cstart, cend = aij->cend,row,col; 101b49de8d1SLois Curfman McInnes int shift = C->indexshift; 102b49de8d1SLois Curfman McInnes 103b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 104b49de8d1SLois Curfman McInnes if (idxm[i] < 0) SETERRQ(1,"MatGetValues_MPIAIJ:Negative row"); 105b49de8d1SLois Curfman McInnes if (idxm[i] >= aij->M) SETERRQ(1,"MatGetValues_MPIAIJ:Row too large"); 106b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 107b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 108b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 109b49de8d1SLois Curfman McInnes if (idxn[j] < 0) SETERRQ(1,"MatGetValues_MPIAIJ:Negative column"); 110b49de8d1SLois Curfman McInnes if (idxn[j] >= aij->N) SETERRQ(1,"MatGetValues_MPIAIJ:Col too large"); 111b49de8d1SLois Curfman McInnes if (idxn[j] >= cstart && idxn[j] < cend){ 112b49de8d1SLois Curfman McInnes col = idxn[j] - cstart; 113b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 114b49de8d1SLois Curfman McInnes } 115b49de8d1SLois Curfman McInnes else { 116b49de8d1SLois Curfman McInnes col = aij->colmap[idxn[j]] + shift; 117b49de8d1SLois Curfman McInnes ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr); 118b49de8d1SLois Curfman McInnes } 119b49de8d1SLois Curfman McInnes } 120b49de8d1SLois Curfman McInnes } 121b49de8d1SLois Curfman McInnes else { 122b49de8d1SLois Curfman McInnes SETERRQ(1,"MatGetValues_MPIAIJ:Only local values currently supported"); 123b49de8d1SLois Curfman McInnes } 124b49de8d1SLois Curfman McInnes } 125b49de8d1SLois Curfman McInnes return 0; 126b49de8d1SLois Curfman McInnes } 127b49de8d1SLois Curfman McInnes 128fc76ce87SLois Curfman McInnes static int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode) 1298a729477SBarry Smith { 13044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 131d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 13217699dbbSLois Curfman McInnes int size = aij->size, *owners = aij->rowners; 13317699dbbSLois Curfman McInnes int rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr; 1341eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 1356abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 1361eb62cbbSBarry Smith InsertMode addv; 1371eb62cbbSBarry Smith Scalar *rvalues,*svalues; 1381eb62cbbSBarry Smith 1391eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 140d65a2f8fSBarry Smith MPI_Allreduce(&aij->insertmode,&addv,1,MPI_INT,MPI_BOR,comm); 141dbb450caSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 142bbb6d6a8SBarry Smith SETERRQ(1,"MatAssemblyBegin_MPIAIJ:Some processors inserted others added"); 1431eb62cbbSBarry Smith } 1441eb62cbbSBarry Smith aij->insertmode = addv; /* in case this processor had no cache */ 1451eb62cbbSBarry Smith 1461eb62cbbSBarry Smith /* first count number of contributors to each processor */ 1470452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 148cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 1490452661fSBarry Smith owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 1501eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1511eb62cbbSBarry Smith idx = aij->stash.idx[i]; 15217699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 1531eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 1541eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 1558a729477SBarry Smith } 1568a729477SBarry Smith } 1578a729477SBarry Smith } 15817699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 1591eb62cbbSBarry Smith 1601eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 1610452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 16217699dbbSLois Curfman McInnes MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm); 16317699dbbSLois Curfman McInnes nreceives = work[rank]; 16417699dbbSLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 16517699dbbSLois Curfman McInnes nmax = work[rank]; 1660452661fSBarry Smith PetscFree(work); 1671eb62cbbSBarry Smith 1681eb62cbbSBarry Smith /* post receives: 1691eb62cbbSBarry Smith 1) each message will consist of ordered pairs 1701eb62cbbSBarry Smith (global index,value) we store the global index as a double 171d6dfbf8fSBarry Smith to simplify the message passing. 1721eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 1731eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 1741eb62cbbSBarry Smith this is a lot of wasted space. 1751eb62cbbSBarry Smith 1761eb62cbbSBarry Smith 1771eb62cbbSBarry Smith This could be done better. 1781eb62cbbSBarry Smith */ 1790452661fSBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 18078b31e54SBarry Smith CHKPTRQ(rvalues); 1810452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request)); 18278b31e54SBarry Smith CHKPTRQ(recv_waits); 1831eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 184416022c9SBarry Smith MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 1851eb62cbbSBarry Smith comm,recv_waits+i); 1861eb62cbbSBarry Smith } 1871eb62cbbSBarry Smith 1881eb62cbbSBarry Smith /* do sends: 1891eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 1901eb62cbbSBarry Smith the ith processor 1911eb62cbbSBarry Smith */ 1920452661fSBarry Smith svalues = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 1930452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 19478b31e54SBarry Smith CHKPTRQ(send_waits); 1950452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 1961eb62cbbSBarry Smith starts[0] = 0; 19717699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1981eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1991eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 2001eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 2011eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 2021eb62cbbSBarry Smith } 2030452661fSBarry Smith PetscFree(owner); 2041eb62cbbSBarry Smith starts[0] = 0; 20517699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2061eb62cbbSBarry Smith count = 0; 20717699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 2081eb62cbbSBarry Smith if (procs[i]) { 209416022c9SBarry Smith MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag, 2101eb62cbbSBarry Smith comm,send_waits+count++); 2111eb62cbbSBarry Smith } 2121eb62cbbSBarry Smith } 2130452661fSBarry Smith PetscFree(starts); PetscFree(nprocs); 2141eb62cbbSBarry Smith 2151eb62cbbSBarry Smith /* Free cache space */ 2162191d07cSBarry Smith PLogInfo(0,"[%d]MatAssemblyBegin_MPIAIJ:Number of off processor values %d\n",rank,aij->stash.n); 21778b31e54SBarry Smith ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr); 2181eb62cbbSBarry Smith 2191eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 2201eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 2211eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 2221eb62cbbSBarry Smith aij->rmax = nmax; 2231eb62cbbSBarry Smith 2241eb62cbbSBarry Smith return 0; 2251eb62cbbSBarry Smith } 22644a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat); 2271eb62cbbSBarry Smith 228fc76ce87SLois Curfman McInnes static int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode) 2291eb62cbbSBarry Smith { 23044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 231dbb450caSBarry Smith Mat_SeqAIJ *C = (Mat_SeqAIJ *) aij->A->data; 2321eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 233416022c9SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr; 234416022c9SBarry Smith int row,col,other_disassembled,shift = C->indexshift; 2351eb62cbbSBarry Smith Scalar *values,val; 2361eb62cbbSBarry Smith InsertMode addv = aij->insertmode; 2371eb62cbbSBarry Smith 2381eb62cbbSBarry Smith /* wait on receives */ 2391eb62cbbSBarry Smith while (count) { 240d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status); 2411eb62cbbSBarry Smith /* unpack receives into our local space */ 242d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 243c60448a5SBarry Smith MPI_Get_count(&recv_status,MPIU_SCALAR,&n); 2441eb62cbbSBarry Smith n = n/3; 2451eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 246227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - aij->rstart; 247227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 2481eb62cbbSBarry Smith val = values[3*i+2]; 2491eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 2501eb62cbbSBarry Smith col -= aij->cstart; 2511eb62cbbSBarry Smith MatSetValues(aij->A,1,&row,1,&col,&val,addv); 2521eb62cbbSBarry Smith } 2531eb62cbbSBarry Smith else { 254227d817aSBarry Smith if (mat->was_assembled) { 255b7c46309SBarry Smith if (!aij->colmap) {ierr = CreateColmap_Private(mat);CHKERRQ(ierr);} 256dbb450caSBarry Smith col = aij->colmap[col] + shift; 257ec8511deSBarry Smith if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) { 2582493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 259227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 260d6dfbf8fSBarry Smith } 2619e25ed09SBarry Smith } 2621eb62cbbSBarry Smith MatSetValues(aij->B,1,&row,1,&col,&val,addv); 2631eb62cbbSBarry Smith } 2641eb62cbbSBarry Smith } 2651eb62cbbSBarry Smith count--; 2661eb62cbbSBarry Smith } 2670452661fSBarry Smith PetscFree(aij->recv_waits); PetscFree(aij->rvalues); 2681eb62cbbSBarry Smith 2691eb62cbbSBarry Smith /* wait on sends */ 2701eb62cbbSBarry Smith if (aij->nsends) { 2710452661fSBarry Smith send_status = (MPI_Status *) PetscMalloc(aij->nsends*sizeof(MPI_Status)); 27278b31e54SBarry Smith CHKPTRQ(send_status); 2731eb62cbbSBarry Smith MPI_Waitall(aij->nsends,aij->send_waits,send_status); 2740452661fSBarry Smith PetscFree(send_status); 2751eb62cbbSBarry Smith } 2760452661fSBarry Smith PetscFree(aij->send_waits); PetscFree(aij->svalues); 2771eb62cbbSBarry Smith 27864119d58SLois Curfman McInnes aij->insertmode = NOT_SET_VALUES; 27978b31e54SBarry Smith ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr); 28078b31e54SBarry Smith ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr); 2811eb62cbbSBarry Smith 2822493cbb0SBarry Smith /* determine if any processor has disassembled, if so we must 2832493cbb0SBarry Smith also disassemble ourselfs, in order that we may reassemble. */ 284227d817aSBarry Smith MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm); 285227d817aSBarry Smith if (mat->was_assembled && !other_disassembled) { 2862493cbb0SBarry Smith ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr); 2872493cbb0SBarry Smith } 2882493cbb0SBarry Smith 289227d817aSBarry Smith if (!mat->was_assembled && mode == FINAL_ASSEMBLY) { 29078b31e54SBarry Smith ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr); 2915e42470aSBarry Smith } 29278b31e54SBarry Smith ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr); 29378b31e54SBarry Smith ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr); 2945e42470aSBarry Smith 2957a0afa10SBarry Smith if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;} 2968a729477SBarry Smith return 0; 2978a729477SBarry Smith } 2988a729477SBarry Smith 2992ee70a88SLois Curfman McInnes static int MatZeroEntries_MPIAIJ(Mat A) 3001eb62cbbSBarry Smith { 30144a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 302dbd7a890SLois Curfman McInnes int ierr; 30378b31e54SBarry Smith ierr = MatZeroEntries(l->A); CHKERRQ(ierr); 30478b31e54SBarry Smith ierr = MatZeroEntries(l->B); CHKERRQ(ierr); 3051eb62cbbSBarry Smith return 0; 3061eb62cbbSBarry Smith } 3071eb62cbbSBarry Smith 3081eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 3091eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 3101eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 3111eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 3121eb62cbbSBarry Smith routine. 3131eb62cbbSBarry Smith */ 31444a69424SLois Curfman McInnes static int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 3151eb62cbbSBarry Smith { 31644a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 31717699dbbSLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 3186abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 31917699dbbSLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 3205392566eSBarry Smith int *rvalues,tag = A->tag,count,base,slen,n,*source; 321d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 322d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 3231eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3241eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 3251eb62cbbSBarry Smith IS istmp; 3261eb62cbbSBarry Smith 32777c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 32878b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 3291eb62cbbSBarry Smith 3301eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3310452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 332cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 3330452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 3341eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3351eb62cbbSBarry Smith idx = rows[i]; 3361eb62cbbSBarry Smith found = 0; 33717699dbbSLois Curfman McInnes for ( j=0; j<size; j++ ) { 3381eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3391eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 3401eb62cbbSBarry Smith } 3411eb62cbbSBarry Smith } 342bbb6d6a8SBarry Smith if (!found) SETERRQ(1,"MatZeroRows_MPIAIJ:Index out of range"); 3431eb62cbbSBarry Smith } 34417699dbbSLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 3451eb62cbbSBarry Smith 3461eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3470452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 34817699dbbSLois Curfman McInnes MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm); 34917699dbbSLois Curfman McInnes nrecvs = work[rank]; 35017699dbbSLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 35117699dbbSLois Curfman McInnes nmax = work[rank]; 3520452661fSBarry Smith PetscFree(work); 3531eb62cbbSBarry Smith 3541eb62cbbSBarry Smith /* post receives: */ 3550452661fSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 35678b31e54SBarry Smith CHKPTRQ(rvalues); 3570452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request)); 35878b31e54SBarry Smith CHKPTRQ(recv_waits); 3591eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 360416022c9SBarry Smith MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i); 3611eb62cbbSBarry Smith } 3621eb62cbbSBarry Smith 3631eb62cbbSBarry Smith /* do sends: 3641eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3651eb62cbbSBarry Smith the ith processor 3661eb62cbbSBarry Smith */ 3670452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 3680452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 36978b31e54SBarry Smith CHKPTRQ(send_waits); 3700452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 3711eb62cbbSBarry Smith starts[0] = 0; 37217699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3731eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3741eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 3751eb62cbbSBarry Smith } 3761eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 3771eb62cbbSBarry Smith 3781eb62cbbSBarry Smith starts[0] = 0; 37917699dbbSLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3801eb62cbbSBarry Smith count = 0; 38117699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 3821eb62cbbSBarry Smith if (procs[i]) { 383416022c9SBarry Smith MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++); 3841eb62cbbSBarry Smith } 3851eb62cbbSBarry Smith } 3860452661fSBarry Smith PetscFree(starts); 3871eb62cbbSBarry Smith 38817699dbbSLois Curfman McInnes base = owners[rank]; 3891eb62cbbSBarry Smith 3901eb62cbbSBarry Smith /* wait on receives */ 3910452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 3921eb62cbbSBarry Smith source = lens + nrecvs; 3931eb62cbbSBarry Smith count = nrecvs; slen = 0; 3941eb62cbbSBarry Smith while (count) { 395d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 3961eb62cbbSBarry Smith /* unpack receives into our local space */ 3971eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_INT,&n); 398d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 399d6dfbf8fSBarry Smith lens[imdex] = n; 4001eb62cbbSBarry Smith slen += n; 4011eb62cbbSBarry Smith count--; 4021eb62cbbSBarry Smith } 4030452661fSBarry Smith PetscFree(recv_waits); 4041eb62cbbSBarry Smith 4051eb62cbbSBarry Smith /* move the data into the send scatter */ 4060452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 4071eb62cbbSBarry Smith count = 0; 4081eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 4091eb62cbbSBarry Smith values = rvalues + i*nmax; 4101eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 4111eb62cbbSBarry Smith lrows[count++] = values[j] - base; 4121eb62cbbSBarry Smith } 4131eb62cbbSBarry Smith } 4140452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 4150452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 4161eb62cbbSBarry Smith 4171eb62cbbSBarry Smith /* actually zap the local rows */ 418416022c9SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 419464493b3SBarry Smith PLogObjectParent(A,istmp); 4200452661fSBarry Smith PetscFree(lrows); 42178b31e54SBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 42278b31e54SBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr); 42378b31e54SBarry Smith ierr = ISDestroy(istmp); CHKERRQ(ierr); 4241eb62cbbSBarry Smith 4251eb62cbbSBarry Smith /* wait on sends */ 4261eb62cbbSBarry Smith if (nsends) { 4270452661fSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status)); 42878b31e54SBarry Smith CHKPTRQ(send_status); 4291eb62cbbSBarry Smith MPI_Waitall(nsends,send_waits,send_status); 4300452661fSBarry Smith PetscFree(send_status); 4311eb62cbbSBarry Smith } 4320452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 4331eb62cbbSBarry Smith 4341eb62cbbSBarry Smith return 0; 4351eb62cbbSBarry Smith } 4361eb62cbbSBarry Smith 437416022c9SBarry Smith static int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy) 4381eb62cbbSBarry Smith { 439416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 4401eb62cbbSBarry Smith int ierr; 441416022c9SBarry Smith 44264119d58SLois Curfman McInnes ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx); CHKERRQ(ierr); 44344cd7ae7SLois Curfman McInnes ierr = MatMult_SeqAIJ(a->A,xx,yy); CHKERRQ(ierr); 44464119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx); CHKERRQ(ierr); 44544cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqAIJ(a->B,a->lvec,yy,yy); CHKERRQ(ierr); 4461eb62cbbSBarry Smith return 0; 4471eb62cbbSBarry Smith } 4481eb62cbbSBarry Smith 449416022c9SBarry Smith static int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 450da3a660dSBarry Smith { 451416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 452da3a660dSBarry Smith int ierr; 45364119d58SLois Curfman McInnes ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr); 45444cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqAIJ(a->A,xx,yy,zz); CHKERRQ(ierr); 45564119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr); 45644cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqAIJ(a->B,a->lvec,zz,zz); CHKERRQ(ierr); 457da3a660dSBarry Smith return 0; 458da3a660dSBarry Smith } 459da3a660dSBarry Smith 460416022c9SBarry Smith static int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy) 461da3a660dSBarry Smith { 462416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 463da3a660dSBarry Smith int ierr; 464da3a660dSBarry Smith 465da3a660dSBarry Smith /* do nondiagonal part */ 46644cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqAIJ(a->B,xx,a->lvec); CHKERRQ(ierr); 467da3a660dSBarry Smith /* send it on its way */ 468416022c9SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES, 46964119d58SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 470da3a660dSBarry Smith /* do local part */ 47144cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqAIJ(a->A,xx,yy); CHKERRQ(ierr); 472da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 473da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 474da3a660dSBarry Smith /* but is not perhaps always true. */ 475416022c9SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES, 47664119d58SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 477da3a660dSBarry Smith return 0; 478da3a660dSBarry Smith } 479da3a660dSBarry Smith 480416022c9SBarry Smith static int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz) 481da3a660dSBarry Smith { 482416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 483da3a660dSBarry Smith int ierr; 484da3a660dSBarry Smith 485da3a660dSBarry Smith /* do nondiagonal part */ 48644cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqAIJ(a->B,xx,a->lvec); CHKERRQ(ierr); 487da3a660dSBarry Smith /* send it on its way */ 488416022c9SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES, 48964119d58SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 490da3a660dSBarry Smith /* do local part */ 49144cd7ae7SLois Curfman McInnes ierr = MatMultTransAdd_SeqAIJ(a->A,xx,yy,zz); CHKERRQ(ierr); 492da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 493da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 494da3a660dSBarry Smith /* but is not perhaps always true. */ 495416022c9SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES, 49664119d58SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 497da3a660dSBarry Smith return 0; 498da3a660dSBarry Smith } 499da3a660dSBarry Smith 5001eb62cbbSBarry Smith /* 5011eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 5021eb62cbbSBarry Smith diagonal block 5031eb62cbbSBarry Smith */ 504416022c9SBarry Smith static int MatGetDiagonal_MPIAIJ(Mat A,Vec v) 5051eb62cbbSBarry Smith { 506416022c9SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 50744cd7ae7SLois Curfman McInnes if (a->M != a->N) 50844cd7ae7SLois Curfman McInnes SETERRQ(1,"MatGetDiagonal_MPIAIJ:Supports only square matrix where A->A is diag block"); 509416022c9SBarry Smith return MatGetDiagonal(a->A,v); 5101eb62cbbSBarry Smith } 5111eb62cbbSBarry Smith 512052efed2SBarry Smith static int MatScale_MPIAIJ(Scalar *aa,Mat A) 513052efed2SBarry Smith { 514052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 515052efed2SBarry Smith int ierr; 516052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 517052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 518052efed2SBarry Smith return 0; 519052efed2SBarry Smith } 520052efed2SBarry Smith 52144a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj) 5221eb62cbbSBarry Smith { 5231eb62cbbSBarry Smith Mat mat = (Mat) obj; 52444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 5251eb62cbbSBarry Smith int ierr; 526a5a9c739SBarry Smith #if defined(PETSC_LOG) 5276e35fa57SLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N); 528a5a9c739SBarry Smith #endif 5290452661fSBarry Smith PetscFree(aij->rowners); 53078b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 53178b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 5320452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 5330452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 5341eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 535a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 5367a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 5370452661fSBarry Smith PetscFree(aij); 538a5a9c739SBarry Smith PLogObjectDestroy(mat); 5390452661fSBarry Smith PetscHeaderDestroy(mat); 5401eb62cbbSBarry Smith return 0; 5411eb62cbbSBarry Smith } 5427857610eSBarry Smith #include "draw.h" 543b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 544ee50ffe9SBarry Smith 545416022c9SBarry Smith static int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 5461eb62cbbSBarry Smith { 547416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 548416022c9SBarry Smith int ierr; 549416022c9SBarry Smith 55017699dbbSLois Curfman McInnes if (aij->size == 1) { 551416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 552416022c9SBarry Smith } 553a523beb4SLois Curfman McInnes else SETERRQ(1,"MatView_MPIAIJ_Binary:Only uniprocessor output supported"); 554416022c9SBarry Smith return 0; 555416022c9SBarry Smith } 556416022c9SBarry Smith 557d7e8b826SBarry Smith static int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer) 558416022c9SBarry Smith { 55944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 560dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 561a56f8943SBarry Smith int ierr, format,shift = C->indexshift,rank; 562d636dbe3SBarry Smith FILE *fd; 56319bcc07fSBarry Smith ViewerType vtype; 564416022c9SBarry Smith 56519bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 56619bcc07fSBarry Smith if (vtype == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) { 56790ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 56890ace30eSBarry Smith if (format == ASCII_FORMAT_INFO_DETAILED) { 56995e01e2fSLois Curfman McInnes int nz, nzalloc, mem, flg; 570a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 57190ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 572a56f8943SBarry Smith ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); 57395e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 57477c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 57595e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 57695e01e2fSLois Curfman McInnes rank,aij->m,nz,nzalloc,mem); 57795e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 57895e01e2fSLois Curfman McInnes rank,aij->m,nz,nzalloc,mem); 57908480c60SBarry Smith ierr = MatGetInfo(aij->A,MAT_LOCAL,&nz,&nzalloc,&mem); 58095e01e2fSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,nz); 58108480c60SBarry Smith ierr = MatGetInfo(aij->B,MAT_LOCAL,&nz,&nzalloc,&mem); 58295e01e2fSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,nz); 583a56f8943SBarry Smith fflush(fd); 58477c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 585a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 586416022c9SBarry Smith return 0; 587416022c9SBarry Smith } 58890ace30eSBarry Smith else if (format == ASCII_FORMAT_INFO) { 58908480c60SBarry Smith return 0; 59008480c60SBarry Smith } 591416022c9SBarry Smith } 592416022c9SBarry Smith 59319bcc07fSBarry Smith if (vtype == DRAW_VIEWER) { 59419bcc07fSBarry Smith Draw draw; 59519bcc07fSBarry Smith PetscTruth isnull; 59619bcc07fSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 59719bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 59819bcc07fSBarry Smith } 59919bcc07fSBarry Smith 60019bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 60190ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 60277c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 603d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 60417699dbbSLois Curfman McInnes aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 6051eb62cbbSBarry Smith aij->cend); 60678b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 60778b31e54SBarry Smith ierr = MatView(aij->B,viewer); CHKERRQ(ierr); 608d13ab20cSBarry Smith fflush(fd); 60977c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 610d13ab20cSBarry Smith } 611416022c9SBarry Smith else { 612a56f8943SBarry Smith int size = aij->size; 613a56f8943SBarry Smith rank = aij->rank; 61417699dbbSLois Curfman McInnes if (size == 1) { 61578b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 61695373324SBarry Smith } 61795373324SBarry Smith else { 61895373324SBarry Smith /* assemble the entire matrix onto first processor. */ 61995373324SBarry Smith Mat A; 620ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 6212eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 62295373324SBarry Smith Scalar *a; 6232ee70a88SLois Curfman McInnes 62417699dbbSLois Curfman McInnes if (!rank) { 625b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 626c451ab25SLois Curfman McInnes CHKERRQ(ierr); 62795373324SBarry Smith } 62895373324SBarry Smith else { 629b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 630c451ab25SLois Curfman McInnes CHKERRQ(ierr); 63195373324SBarry Smith } 632464493b3SBarry Smith PLogObjectParent(mat,A); 633416022c9SBarry Smith 63495373324SBarry Smith /* copy over the A part */ 635ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 6362ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 63795373324SBarry Smith row = aij->rstart; 638dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 63995373324SBarry Smith for ( i=0; i<m; i++ ) { 640416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 64195373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 64295373324SBarry Smith } 6432ee70a88SLois Curfman McInnes aj = Aloc->j; 644dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 64595373324SBarry Smith 64695373324SBarry Smith /* copy over the B part */ 647ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 6482ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 64995373324SBarry Smith row = aij->rstart; 6500452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 651dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 65295373324SBarry Smith for ( i=0; i<m; i++ ) { 653416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 65495373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 65595373324SBarry Smith } 6560452661fSBarry Smith PetscFree(ct); 65778b31e54SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 65878b31e54SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 65917699dbbSLois Curfman McInnes if (!rank) { 66078b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 66195373324SBarry Smith } 66278b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 66395373324SBarry Smith } 66495373324SBarry Smith } 6651eb62cbbSBarry Smith return 0; 6661eb62cbbSBarry Smith } 6671eb62cbbSBarry Smith 668416022c9SBarry Smith static int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 669416022c9SBarry Smith { 670416022c9SBarry Smith Mat mat = (Mat) obj; 671416022c9SBarry Smith int ierr; 67219bcc07fSBarry Smith ViewerType vtype; 673416022c9SBarry Smith 674416022c9SBarry Smith if (!viewer) { 67519bcc07fSBarry Smith viewer = STDOUT_VIEWER_SELF; 676416022c9SBarry Smith } 67719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 67819bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER || 67919bcc07fSBarry Smith vtype == DRAW_VIEWER || vtype == MATLAB_VIEWER) { 680d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 681416022c9SBarry Smith } 68219bcc07fSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 683416022c9SBarry Smith return MatView_MPIAIJ_Binary(mat,viewer); 684416022c9SBarry Smith } 685416022c9SBarry Smith return 0; 686416022c9SBarry Smith } 687416022c9SBarry Smith 6881eb62cbbSBarry Smith /* 6891eb62cbbSBarry Smith This has to provide several versions. 6901eb62cbbSBarry Smith 6911eb62cbbSBarry Smith 1) per sequential 6921eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 6931eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 694d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 6951eb62cbbSBarry Smith */ 696fc76ce87SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 697dbb450caSBarry Smith double fshift,int its,Vec xx) 6988a729477SBarry Smith { 69944a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 700d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 701ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 7026abc6512SBarry Smith Scalar zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts; 7036abc6512SBarry Smith int ierr,*idx, *diag; 704416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 705da3a660dSBarry Smith Vec tt; 7068a729477SBarry Smith 707d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 708dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 709dbb450caSBarry Smith ls = ls + shift; 710ec8511deSBarry Smith if (!A->diag) {if ((ierr = MatMarkDiag_SeqAIJ(AA))) return ierr;} 711d6dfbf8fSBarry Smith diag = A->diag; 712acb40c82SBarry Smith if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) { 71348d91487SBarry Smith SETERRQ(1,"MatRelax_MPIAIJ:Option not supported"); 714acb40c82SBarry Smith } 715da3a660dSBarry Smith if (flag & SOR_EISENSTAT) { 716da3a660dSBarry Smith /* Let A = L + U + D; where L is lower trianglar, 717da3a660dSBarry Smith U is upper triangular, E is diagonal; This routine applies 718da3a660dSBarry Smith 719da3a660dSBarry Smith (L + E)^{-1} A (U + E)^{-1} 720da3a660dSBarry Smith 721da3a660dSBarry Smith to a vector efficiently using Eisenstat's trick. This is for 722da3a660dSBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 723da3a660dSBarry Smith is the relaxation factor. 724da3a660dSBarry Smith */ 72578b31e54SBarry Smith ierr = VecDuplicate(xx,&tt); CHKERRQ(ierr); 726464493b3SBarry Smith PLogObjectParent(matin,tt); 727da3a660dSBarry Smith VecGetArray(tt,&t); 728da3a660dSBarry Smith scale = (2.0/omega) - 1.0; 729da3a660dSBarry Smith /* x = (E + U)^{-1} b */ 730da3a660dSBarry Smith VecSet(&zero,mat->lvec); 73164119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 73278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 733da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 734da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 735dbb450caSBarry Smith idx = A->j + diag[i] + !shift; 736dbb450caSBarry Smith v = A->a + diag[i] + !shift; 737da3a660dSBarry Smith sum = b[i]; 738da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 739dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 740da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 741dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 742dbb450caSBarry Smith v = B->a + B->i[i] + shift; 743da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 744da3a660dSBarry Smith x[i] = omega*(sum/d); 745da3a660dSBarry Smith } 74664119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 74778b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 748da3a660dSBarry Smith 749da3a660dSBarry Smith /* t = b - (2*E - D)x */ 750da3a660dSBarry Smith v = A->a; 751dbb450caSBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 752da3a660dSBarry Smith 753da3a660dSBarry Smith /* t = (E + L)^{-1}t */ 754dbb450caSBarry Smith ts = t + shift; /* shifted by one for index start of a or mat->j*/ 755da3a660dSBarry Smith diag = A->diag; 756da3a660dSBarry Smith VecSet(&zero,mat->lvec); 75764119d58SLois Curfman McInnes ierr = VecPipelineBegin(tt,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 75878b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 759da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 760da3a660dSBarry Smith n = diag[i] - A->i[i]; 761dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 762dbb450caSBarry Smith v = A->a + A->i[i] + shift; 763da3a660dSBarry Smith sum = t[i]; 764da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 765dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 766da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 767dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 768dbb450caSBarry Smith v = B->a + B->i[i] + shift; 769da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 770da3a660dSBarry Smith t[i] = omega*(sum/d); 771da3a660dSBarry Smith } 77264119d58SLois Curfman McInnes ierr = VecPipelineEnd(tt,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 77378b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 774da3a660dSBarry Smith /* x = x + t */ 775da3a660dSBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 776da3a660dSBarry Smith VecDestroy(tt); 777da3a660dSBarry Smith return 0; 778da3a660dSBarry Smith } 779da3a660dSBarry Smith 7801eb62cbbSBarry Smith 781d6dfbf8fSBarry Smith if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){ 782da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 783da3a660dSBarry Smith VecSet(&zero,mat->lvec); VecSet(&zero,xx); 784da3a660dSBarry Smith } 785da3a660dSBarry Smith else { 78664119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 78778b31e54SBarry Smith CHKERRQ(ierr); 78864119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 78978b31e54SBarry Smith CHKERRQ(ierr); 790da3a660dSBarry Smith } 791d6dfbf8fSBarry Smith while (its--) { 792d6dfbf8fSBarry Smith /* go down through the rows */ 79364119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 79478b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 795d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 796d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 797dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 798dbb450caSBarry Smith v = A->a + A->i[i] + shift; 799d6dfbf8fSBarry Smith sum = b[i]; 800d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 801dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 802d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 803dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 804dbb450caSBarry Smith v = B->a + B->i[i] + shift; 805d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 806d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 807d6dfbf8fSBarry Smith } 80864119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 80978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 810d6dfbf8fSBarry Smith /* come up through the rows */ 81164119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 81278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 813d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 814d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 815dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 816dbb450caSBarry Smith v = A->a + A->i[i] + shift; 817d6dfbf8fSBarry Smith sum = b[i]; 818d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 819dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 820d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 821dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 822dbb450caSBarry Smith v = B->a + B->i[i] + shift; 823d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 824d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 825d6dfbf8fSBarry Smith } 82664119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 82778b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 828d6dfbf8fSBarry Smith } 829d6dfbf8fSBarry Smith } 830d6dfbf8fSBarry Smith else if (flag & SOR_FORWARD_SWEEP){ 831da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 832da3a660dSBarry Smith VecSet(&zero,mat->lvec); 83364119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 83478b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 835da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 836da3a660dSBarry Smith n = diag[i] - A->i[i]; 837dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 838dbb450caSBarry Smith v = A->a + A->i[i] + shift; 839da3a660dSBarry Smith sum = b[i]; 840da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 841dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 842da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 843dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 844dbb450caSBarry Smith v = B->a + B->i[i] + shift; 845da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 846da3a660dSBarry Smith x[i] = omega*(sum/d); 847da3a660dSBarry Smith } 84864119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 84978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 850da3a660dSBarry Smith its--; 851da3a660dSBarry Smith } 852d6dfbf8fSBarry Smith while (its--) { 85364119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 85478b31e54SBarry Smith CHKERRQ(ierr); 85564119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 85678b31e54SBarry Smith CHKERRQ(ierr); 85764119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 85878b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 859d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 860d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 861dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 862dbb450caSBarry Smith v = A->a + A->i[i] + shift; 863d6dfbf8fSBarry Smith sum = b[i]; 864d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 865dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 866d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 867dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 868dbb450caSBarry Smith v = B->a + B->i[i] + shift; 869d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 870d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 871d6dfbf8fSBarry Smith } 87264119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 87378b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 874d6dfbf8fSBarry Smith } 875d6dfbf8fSBarry Smith } 876d6dfbf8fSBarry Smith else if (flag & SOR_BACKWARD_SWEEP){ 877da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 878da3a660dSBarry Smith VecSet(&zero,mat->lvec); 87964119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 88078b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 881da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 882da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 883dbb450caSBarry Smith idx = A->j + diag[i] + !shift; 884dbb450caSBarry Smith v = A->a + diag[i] + !shift; 885da3a660dSBarry Smith sum = b[i]; 886da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 887dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 888da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 889dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 890dbb450caSBarry Smith v = B->a + B->i[i] + shift; 891da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 892da3a660dSBarry Smith x[i] = omega*(sum/d); 893da3a660dSBarry Smith } 89464119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 89578b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 896da3a660dSBarry Smith its--; 897da3a660dSBarry Smith } 898d6dfbf8fSBarry Smith while (its--) { 89964119d58SLois Curfman McInnes ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_DOWN, 90078b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 90164119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_DOWN, 90278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 90364119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 90478b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 905d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 906d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 907dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 908dbb450caSBarry Smith v = A->a + A->i[i] + shift; 909d6dfbf8fSBarry Smith sum = b[i]; 910d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 911dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 912d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 913dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 914dbb450caSBarry Smith v = B->a + B->i[i] + shift; 915d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 916d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 917d6dfbf8fSBarry Smith } 91864119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 91978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 920d6dfbf8fSBarry Smith } 921d6dfbf8fSBarry Smith } 922d6dfbf8fSBarry Smith else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 923da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 92444cd7ae7SLois Curfman McInnes return MatRelax_SeqAIJ(mat->A,bb,omega,flag,fshift,its,xx); 925da3a660dSBarry Smith } 92664119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 92778b31e54SBarry Smith CHKERRQ(ierr); 92864119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 92978b31e54SBarry Smith CHKERRQ(ierr); 930d6dfbf8fSBarry Smith while (its--) { 931d6dfbf8fSBarry Smith /* go down through the rows */ 932d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 933d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 934dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 935dbb450caSBarry Smith v = A->a + A->i[i] + shift; 936d6dfbf8fSBarry Smith sum = b[i]; 937d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 938dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 939d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 940dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 941dbb450caSBarry Smith v = B->a + B->i[i] + shift; 942d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 943d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 944d6dfbf8fSBarry Smith } 945d6dfbf8fSBarry Smith /* come up through the rows */ 946d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 947d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 948dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 949dbb450caSBarry Smith v = A->a + A->i[i] + shift; 950d6dfbf8fSBarry Smith sum = b[i]; 951d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 952dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 953d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 954dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 955dbb450caSBarry Smith v = B->a + B->i[i] + shift; 956d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 957d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 958d6dfbf8fSBarry Smith } 959d6dfbf8fSBarry Smith } 960d6dfbf8fSBarry Smith } 961d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 962da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 96344cd7ae7SLois Curfman McInnes return MatRelax_SeqAIJ(mat->A,bb,omega,flag,fshift,its,xx); 964da3a660dSBarry Smith } 96564119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 96678b31e54SBarry Smith CHKERRQ(ierr); 96764119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 96878b31e54SBarry Smith CHKERRQ(ierr); 969d6dfbf8fSBarry Smith while (its--) { 970d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 971d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 972dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 973dbb450caSBarry Smith v = A->a + A->i[i] + shift; 974d6dfbf8fSBarry Smith sum = b[i]; 975d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 976dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 977d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 978dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 979dbb450caSBarry Smith v = B->a + B->i[i] + shift; 980d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 981d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 982d6dfbf8fSBarry Smith } 983d6dfbf8fSBarry Smith } 984d6dfbf8fSBarry Smith } 985d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 986da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 98744cd7ae7SLois Curfman McInnes return MatRelax_SeqAIJ(mat->A,bb,omega,flag,fshift,its,xx); 988da3a660dSBarry Smith } 98964119d58SLois Curfman McInnes ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL, 99078b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 99164119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL, 99278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 993d6dfbf8fSBarry Smith while (its--) { 994d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 995d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 996dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 997dbb450caSBarry Smith v = A->a + A->i[i] + shift; 998d6dfbf8fSBarry Smith sum = b[i]; 999d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 1000dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 1001d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 1002dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 1003dbb450caSBarry Smith v = B->a + B->i[i] + shift; 1004d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 1005d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 1006d6dfbf8fSBarry Smith } 1007d6dfbf8fSBarry Smith } 1008d6dfbf8fSBarry Smith } 10098a729477SBarry Smith return 0; 10108a729477SBarry Smith } 1011a66be287SLois Curfman McInnes 1012fc76ce87SLois Curfman McInnes static int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,int *nz, 1013fc76ce87SLois Curfman McInnes int *nzalloc,int *mem) 1014a66be287SLois Curfman McInnes { 1015a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1016a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 1017a66be287SLois Curfman McInnes int ierr, isend[3], irecv[3], nzA, nzallocA, memA; 1018a66be287SLois Curfman McInnes 101978b31e54SBarry Smith ierr = MatGetInfo(A,MAT_LOCAL,&nzA,&nzallocA,&memA); CHKERRQ(ierr); 102078b31e54SBarry Smith ierr = MatGetInfo(B,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr); 1021a66be287SLois Curfman McInnes isend[0] += nzA; isend[1] += nzallocA; isend[2] += memA; 1022a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 1023bcd2baecSBarry Smith if (nz) *nz = isend[0]; 1024bcd2baecSBarry Smith if (nzalloc) *nzalloc = isend[1]; 1025bcd2baecSBarry Smith if (mem) *mem = isend[2]; 1026a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 1027d65a2f8fSBarry Smith MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_MAX,matin->comm); 1028bcd2baecSBarry Smith if (nz) *nz = irecv[0]; 1029bcd2baecSBarry Smith if (nzalloc) *nzalloc = irecv[1]; 1030bcd2baecSBarry Smith if (mem) *mem = irecv[2]; 1031a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 1032d65a2f8fSBarry Smith MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_SUM,matin->comm); 1033bcd2baecSBarry Smith if (nz) *nz = irecv[0]; 1034bcd2baecSBarry Smith if (nzalloc) *nzalloc = irecv[1]; 1035bcd2baecSBarry Smith if (mem) *mem = irecv[2]; 1036a66be287SLois Curfman McInnes } 1037a66be287SLois Curfman McInnes return 0; 1038a66be287SLois Curfman McInnes } 1039a66be287SLois Curfman McInnes 1040299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*); 1041299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*); 1042299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double); 1043299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *); 1044299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec); 1045299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec); 1046299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec); 1047299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec); 1048299609e3SLois Curfman McInnes 1049416022c9SBarry Smith static int MatSetOption_MPIAIJ(Mat A,MatOption op) 1050c74985f6SBarry Smith { 1051c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1052c74985f6SBarry Smith 1053c0bbcb79SLois Curfman McInnes if (op == NO_NEW_NONZERO_LOCATIONS || 1054c0bbcb79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 1055c0bbcb79SLois Curfman McInnes op == COLUMNS_SORTED || 1056c0bbcb79SLois Curfman McInnes op == ROW_ORIENTED) { 1057c0bbcb79SLois Curfman McInnes MatSetOption(a->A,op); 1058c0bbcb79SLois Curfman McInnes MatSetOption(a->B,op); 1059c74985f6SBarry Smith } 1060c0bbcb79SLois Curfman McInnes else if (op == ROWS_SORTED || 1061c0bbcb79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 1062c0bbcb79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 1063c0bbcb79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 106494a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_MPIAIJ:Option ignored\n"); 10654b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) { 10664b0e389bSBarry Smith a->roworiented = 0; 10674b0e389bSBarry Smith MatSetOption(a->A,op); 10684b0e389bSBarry Smith MatSetOption(a->B,op); 10694b0e389bSBarry Smith } 1070c0bbcb79SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 10714d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:NO_NEW_DIAGONALS");} 1072c0bbcb79SLois Curfman McInnes else 10734d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:unknown option");} 1074c74985f6SBarry Smith return 0; 1075c74985f6SBarry Smith } 1076c74985f6SBarry Smith 1077d1710a03SLois Curfman McInnes static int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1078c74985f6SBarry Smith { 107944a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1080c74985f6SBarry Smith *m = mat->M; *n = mat->N; 1081c74985f6SBarry Smith return 0; 1082c74985f6SBarry Smith } 1083c74985f6SBarry Smith 1084d1710a03SLois Curfman McInnes static int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1085c74985f6SBarry Smith { 108644a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1087b7c46309SBarry Smith *m = mat->m; *n = mat->N; 1088c74985f6SBarry Smith return 0; 1089c74985f6SBarry Smith } 1090c74985f6SBarry Smith 1091d1710a03SLois Curfman McInnes static int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1092c74985f6SBarry Smith { 109344a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1094c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 1095c74985f6SBarry Smith return 0; 1096c74985f6SBarry Smith } 1097c74985f6SBarry Smith 10986d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 10996d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 11006d84be18SBarry Smith 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 11097a0afa10SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,"MatGetRow_MPIAIJ:Already active"); 11107a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11117a0afa10SBarry Smith 111270f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11137a0afa10SBarry Smith /* 11147a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11157a0afa10SBarry Smith */ 11167a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 11177a0afa10SBarry Smith int max = 1,n = mat->n,tmp; 11187a0afa10SBarry Smith for ( i=0; i<n; i++ ) { 11197a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11207a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11217a0afa10SBarry Smith } 11227a0afa10SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar))); 11237a0afa10SBarry Smith CHKPTRQ(mat->rowvalues); 11247a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11257a0afa10SBarry Smith } 11267a0afa10SBarry Smith 11277a0afa10SBarry Smith 1128b49de8d1SLois Curfman McInnes if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIAIJ:Only local rows") 1129abc0e9e4SLois Curfman McInnes lrow = row - rstart; 113039e00950SLois Curfman McInnes 1131154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1132154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1133154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 11346d84be18SBarry Smith ierr = MatGetRow_SeqAIJ(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 11356d84be18SBarry Smith ierr = MatGetRow_SeqAIJ(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1136154123eaSLois Curfman McInnes nztot = nzA + nzB; 1137154123eaSLois Curfman McInnes 113870f0671dSBarry Smith cmap = mat->garray; 1139154123eaSLois Curfman McInnes if (v || idx) { 1140154123eaSLois Curfman McInnes if (nztot) { 1141154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 114270f0671dSBarry Smith int imark = -1; 1143154123eaSLois Curfman McInnes if (v) { 114470f0671dSBarry Smith *v = v_p = mat->rowvalues; 114539e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 114670f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1147154123eaSLois Curfman McInnes else break; 1148154123eaSLois Curfman McInnes } 1149154123eaSLois Curfman McInnes imark = i; 115070f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 115170f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1152154123eaSLois Curfman McInnes } 1153154123eaSLois Curfman McInnes if (idx) { 115470f0671dSBarry Smith *idx = idx_p = mat->rowindices; 115570f0671dSBarry Smith if (imark > -1) { 115670f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 115770f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 115870f0671dSBarry Smith } 115970f0671dSBarry Smith } else { 1160154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 116170f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1162154123eaSLois Curfman McInnes else break; 1163154123eaSLois Curfman McInnes } 1164154123eaSLois Curfman McInnes imark = i; 116570f0671dSBarry Smith } 116670f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 116770f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 116839e00950SLois Curfman McInnes } 116939e00950SLois Curfman McInnes } 117039e00950SLois Curfman McInnes else {*idx = 0; *v=0;} 1171154123eaSLois Curfman McInnes } 117239e00950SLois Curfman McInnes *nz = nztot; 11736d84be18SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 11746d84be18SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 117539e00950SLois Curfman McInnes return 0; 117639e00950SLois Curfman McInnes } 117739e00950SLois Curfman McInnes 11786d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 117939e00950SLois Curfman McInnes { 11807a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 11817a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 11827a0afa10SBarry Smith SETERRQ(1,"MatRestoreRow_MPIAIJ:MatGetRow not called"); 11837a0afa10SBarry Smith } 11847a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 118539e00950SLois Curfman McInnes return 0; 118639e00950SLois Curfman McInnes } 118739e00950SLois Curfman McInnes 1188cddf8d76SBarry Smith static int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1189855ac2c5SLois Curfman McInnes { 1190855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1191ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1192416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1193416022c9SBarry Smith double sum = 0.0; 119404ca555eSLois Curfman McInnes Scalar *v; 119504ca555eSLois Curfman McInnes 119617699dbbSLois Curfman McInnes if (aij->size == 1) { 119714183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 119837fa93a5SLois Curfman McInnes } else { 119904ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 120004ca555eSLois Curfman McInnes v = amat->a; 120104ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 120204ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 120304ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 120404ca555eSLois Curfman McInnes #else 120504ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 120604ca555eSLois Curfman McInnes #endif 120704ca555eSLois Curfman McInnes } 120804ca555eSLois Curfman McInnes v = bmat->a; 120904ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 121004ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 121104ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 121204ca555eSLois Curfman McInnes #else 121304ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 121404ca555eSLois Curfman McInnes #endif 121504ca555eSLois Curfman McInnes } 12166d84be18SBarry Smith MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm); 121704ca555eSLois Curfman McInnes *norm = sqrt(*norm); 121804ca555eSLois Curfman McInnes } 121904ca555eSLois Curfman McInnes else if (type == NORM_1) { /* max column norm */ 122004ca555eSLois Curfman McInnes double *tmp, *tmp2; 122104ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 12220452661fSBarry Smith tmp = (double *) PetscMalloc( aij->N*sizeof(double) ); CHKPTRQ(tmp); 12230452661fSBarry Smith tmp2 = (double *) PetscMalloc( aij->N*sizeof(double) ); CHKPTRQ(tmp2); 1224cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 122504ca555eSLois Curfman McInnes *norm = 0.0; 122604ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 122704ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1228579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 122904ca555eSLois Curfman McInnes } 123004ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 123104ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1232579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 123304ca555eSLois Curfman McInnes } 12346d84be18SBarry Smith MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm); 123504ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 123604ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 123704ca555eSLois Curfman McInnes } 12380452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 123904ca555eSLois Curfman McInnes } 124004ca555eSLois Curfman McInnes else if (type == NORM_INFINITY) { /* max row norm */ 1241515d9167SLois Curfman McInnes double ntemp = 0.0; 124204ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1243dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 124404ca555eSLois Curfman McInnes sum = 0.0; 124504ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1246cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 124704ca555eSLois Curfman McInnes } 1248dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 124904ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1250cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 125104ca555eSLois Curfman McInnes } 1252515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 125304ca555eSLois Curfman McInnes } 12546d84be18SBarry Smith MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm); 125504ca555eSLois Curfman McInnes } 125604ca555eSLois Curfman McInnes else { 1257515d9167SLois Curfman McInnes SETERRQ(1,"MatNorm_MPIAIJ:No support for two norm"); 125804ca555eSLois Curfman McInnes } 125937fa93a5SLois Curfman McInnes } 1260855ac2c5SLois Curfman McInnes return 0; 1261855ac2c5SLois Curfman McInnes } 1262855ac2c5SLois Curfman McInnes 12630de55854SLois Curfman McInnes static int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1264b7c46309SBarry Smith { 1265b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1266dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1267416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1268416022c9SBarry Smith Mat B; 1269b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 1270b7c46309SBarry Smith Scalar *array; 1271b7c46309SBarry Smith 12723638b69dSLois Curfman McInnes if (matout == PETSC_NULL && M != N) 12733638b69dSLois Curfman McInnes SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place"); 1274b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,PETSC_NULL,0, 1275b4fd4287SBarry Smith PETSC_NULL,&B); CHKERRQ(ierr); 1276b7c46309SBarry Smith 1277b7c46309SBarry Smith /* copy over the A part */ 1278ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1279b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1280b7c46309SBarry Smith row = a->rstart; 1281dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1282b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1283416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1284b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1285b7c46309SBarry Smith } 1286b7c46309SBarry Smith aj = Aloc->j; 12874af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1288b7c46309SBarry Smith 1289b7c46309SBarry Smith /* copy over the B part */ 1290ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1291b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1292b7c46309SBarry Smith row = a->rstart; 12930452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1294dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1295b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1296416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1297b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1298b7c46309SBarry Smith } 12990452661fSBarry Smith PetscFree(ct); 1300b7c46309SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1301b7c46309SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 13023638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13030de55854SLois Curfman McInnes *matout = B; 13040de55854SLois Curfman McInnes } else { 13050de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13060452661fSBarry Smith PetscFree(a->rowners); 13070de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13080de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 13090452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 13100452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13110de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1312a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13130452661fSBarry Smith PetscFree(a); 1314416022c9SBarry Smith PetscMemcpy(A,B,sizeof(struct _Mat)); 13150452661fSBarry Smith PetscHeaderDestroy(B); 13160de55854SLois Curfman McInnes } 1317b7c46309SBarry Smith return 0; 1318b7c46309SBarry Smith } 1319b7c46309SBarry Smith 1320682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat); 1321682d7d0cSBarry Smith static int MatPrintHelp_MPIAIJ(Mat A) 1322682d7d0cSBarry Smith { 1323682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1324682d7d0cSBarry Smith 1325682d7d0cSBarry Smith if (!a->rank) return MatPrintHelp_SeqAIJ(a->A); 1326682d7d0cSBarry Smith else return 0; 1327682d7d0cSBarry Smith } 1328682d7d0cSBarry Smith 1329fc76ce87SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MatType,Mat *); 13303d1612f7SBarry Smith static int MatConvertSameType_MPIAIJ(Mat,Mat *,int); 13312f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 1332598137ffSSatish Balay int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **); 13338a729477SBarry Smith /* -------------------------------------------------------------------*/ 13342ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ, 133539e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 133644a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 133744a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 1338299609e3SLois Curfman McInnes MatSolve_MPIAIJ,MatSolveAdd_MPIAIJ, 1339299609e3SLois Curfman McInnes MatSolveTrans_MPIAIJ,MatSolveTransAdd_MPIAIJ, 1340299609e3SLois Curfman McInnes MatLUFactor_MPIAIJ,0, 134144a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1342b7c46309SBarry Smith MatTranspose_MPIAIJ, 1343a66be287SLois Curfman McInnes MatGetInfo_MPIAIJ,0, 1344855ac2c5SLois Curfman McInnes MatGetDiagonal_MPIAIJ,0,MatNorm_MPIAIJ, 1345ee50ffe9SBarry Smith MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ, 13461eb62cbbSBarry Smith 0, 1347299609e3SLois Curfman McInnes MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ, 1348299609e3SLois Curfman McInnes MatGetReordering_MPIAIJ, 1349299609e3SLois Curfman McInnes MatLUFactorSymbolic_MPIAIJ,MatLUFactorNumeric_MPIAIJ,0,0, 1350d1710a03SLois Curfman McInnes MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ, 1351299609e3SLois Curfman McInnes MatILUFactorSymbolic_MPIAIJ,0, 13523d1612f7SBarry Smith 0,0,MatConvert_MPIAIJ,0,0,MatConvertSameType_MPIAIJ,0,0, 1353b49de8d1SLois Curfman McInnes 0,0,0, 1354598137ffSSatish Balay MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0, 1355052efed2SBarry Smith MatPrintHelp_MPIAIJ, 1356052efed2SBarry Smith MatScale_MPIAIJ}; 13578a729477SBarry Smith 13581987afe7SBarry Smith /*@C 1359ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 13603a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 13613a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 13623a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 13633a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 13648a729477SBarry Smith 13658a729477SBarry Smith Input Parameters: 13661eb62cbbSBarry Smith . comm - MPI communicator 13677d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 13687d3e4905SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated 13697d3e4905SLois Curfman McInnes if N is given) 13707d3e4905SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 13717d3e4905SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated 13727d3e4905SLois Curfman McInnes if n is given) 1373ab693e5aSLois Curfman McInnes . d_nz - number of nonzeros per row in diagonal portion of local submatrix 1374ff756334SLois Curfman McInnes (same for all local rows) 1375ab693e5aSLois Curfman McInnes . d_nzz - number of nonzeros per row in diagonal portion of local submatrix 1376ab693e5aSLois Curfman McInnes or null (possibly different for each row). You must leave room 1377ab693e5aSLois Curfman McInnes for the diagonal entry even if it is zero. 1378ab693e5aSLois Curfman McInnes . o_nz - number of nonzeros per row in off-diagonal portion of local 1379ab693e5aSLois Curfman McInnes submatrix (same for all local rows). 1380ab693e5aSLois Curfman McInnes . o_nzz - number of nonzeros per row in off-diagonal portion of local 1381ab693e5aSLois Curfman McInnes submatrix or null (possibly different for each row). 13828a729477SBarry Smith 1383ff756334SLois Curfman McInnes Output Parameter: 138444cd7ae7SLois Curfman McInnes . A - the matrix 13858a729477SBarry Smith 1386ff756334SLois Curfman McInnes Notes: 1387ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1388ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 13890002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 13900002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1391ff756334SLois Curfman McInnes 1392ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1393ff756334SLois Curfman McInnes (possibly both). 1394ff756334SLois Curfman McInnes 13955511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 13965511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 13975511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 13985511cfe3SLois Curfman McInnes 13995511cfe3SLois Curfman McInnes Options Database Keys: 14005511cfe3SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 14016e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 14026e62573dSLois Curfman McInnes $ (max limit=5) 14036e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 14046e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 14056e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 14065511cfe3SLois Curfman McInnes 1407e0245417SLois Curfman McInnes Storage Information: 1408e0245417SLois Curfman McInnes For a square global matrix we define each processor's diagonal portion 1409e0245417SLois Curfman McInnes to be its local rows and the corresponding columns (a square submatrix); 1410e0245417SLois Curfman McInnes each processor's off-diagonal portion encompasses the remainder of the 1411e0245417SLois Curfman McInnes local matrix (a rectangular submatrix). 1412e0245417SLois Curfman McInnes 1413e0245417SLois Curfman McInnes The user can specify preallocated storage for the diagonal part of 14145ace5be8SLois Curfman McInnes the local submatrix with either d_nz or d_nnz (not both). Set 14155ace5be8SLois Curfman McInnes d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic 14165ace5be8SLois Curfman McInnes memory allocation. Likewise, specify preallocated storage for the 14175ace5be8SLois Curfman McInnes off-diagonal part of the local submatrix with o_nz or o_nnz (not both). 1418ff756334SLois Curfman McInnes 14195511cfe3SLois Curfman McInnes Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In 14205511cfe3SLois Curfman McInnes the figure below we depict these three local rows and all columns (0-11). 14212191d07cSBarry Smith 1422b810aeb4SBarry Smith $ 0 1 2 3 4 5 6 7 8 9 10 11 1423b810aeb4SBarry Smith $ ------------------- 14245511cfe3SLois Curfman McInnes $ row 3 | o o o d d d o o o o o o 14255511cfe3SLois Curfman McInnes $ row 4 | o o o d d d o o o o o o 14265511cfe3SLois Curfman McInnes $ row 5 | o o o d d d o o o o o o 14275511cfe3SLois Curfman McInnes $ ------------------- 1428b810aeb4SBarry Smith $ 1429b810aeb4SBarry Smith 14305511cfe3SLois Curfman McInnes Thus, any entries in the d locations are stored in the d (diagonal) 14315511cfe3SLois Curfman McInnes submatrix, and any entries in the o locations are stored in the 14325511cfe3SLois Curfman McInnes o (off-diagonal) submatrix. Note that the d and the o submatrices are 14335511cfe3SLois Curfman McInnes stored simply in the MATSEQAIJ format for compressed row storage. 14345511cfe3SLois Curfman McInnes 14355511cfe3SLois Curfman McInnes Now d_nz should indicate the number of nonzeros per row in the d matrix, 14365511cfe3SLois Curfman McInnes and o_nz should indicate the number of nonzeros per row in the o matrix. 14375511cfe3SLois Curfman McInnes In general, for PDE problems in which most nonzeros are near the diagonal, 14383a511b96SLois Curfman McInnes one expects d_nz >> o_nz. For additional details, see the users manual 14393a511b96SLois Curfman McInnes chapter on matrices and the file $(PETSC_DIR)/Performance. 14403a511b96SLois Curfman McInnes 1441dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1442ff756334SLois Curfman McInnes 1443fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 14448a729477SBarry Smith @*/ 14451eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 144644cd7ae7SLois Curfman McInnes int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A) 14478a729477SBarry Smith { 144844cd7ae7SLois Curfman McInnes Mat B; 144944cd7ae7SLois Curfman McInnes Mat_MPIAIJ *b; 14506abc6512SBarry Smith int ierr, i,sum[2],work[2]; 1451416022c9SBarry Smith 145244cd7ae7SLois Curfman McInnes *A = 0; 145344cd7ae7SLois Curfman McInnes PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATMPIAIJ,comm); 145444cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 145544cd7ae7SLois Curfman McInnes B->data = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b); 145644cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_MPIAIJ)); 145744cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 145844cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_MPIAIJ; 145944cd7ae7SLois Curfman McInnes B->view = MatView_MPIAIJ; 146044cd7ae7SLois Curfman McInnes B->factor = 0; 146144cd7ae7SLois Curfman McInnes B->assembled = PETSC_FALSE; 1462d6dfbf8fSBarry Smith 146344cd7ae7SLois Curfman McInnes b->insertmode = NOT_SET_VALUES; 146444cd7ae7SLois Curfman McInnes MPI_Comm_rank(comm,&b->rank); 146544cd7ae7SLois Curfman McInnes MPI_Comm_size(comm,&b->size); 14661eb62cbbSBarry Smith 1467b4fd4287SBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) 14682e0155e0SLois Curfman McInnes SETERRQ(1,"MatCreateMPIAIJ:Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz"); 14691987afe7SBarry Smith 1470dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE || N == PETSC_DECIDE) { 14711eb62cbbSBarry Smith work[0] = m; work[1] = n; 1472d65a2f8fSBarry Smith MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm ); 1473dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE) M = sum[0]; 1474dbd7a890SLois Curfman McInnes if (N == PETSC_DECIDE) N = sum[1]; 14751eb62cbbSBarry Smith } 147644cd7ae7SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);} 147744cd7ae7SLois Curfman McInnes if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);} 147844cd7ae7SLois Curfman McInnes b->m = m; B->m = m; 147944cd7ae7SLois Curfman McInnes b->n = n; B->n = n; 148044cd7ae7SLois Curfman McInnes b->N = N; B->N = N; 148144cd7ae7SLois Curfman McInnes b->M = M; B->M = M; 14821eb62cbbSBarry Smith 14831eb62cbbSBarry Smith /* build local table of row and column ownerships */ 148444cd7ae7SLois Curfman McInnes b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners); 148544cd7ae7SLois Curfman McInnes PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ)); 148644cd7ae7SLois Curfman McInnes b->cowners = b->rowners + b->size + 1; 148744cd7ae7SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm); 148844cd7ae7SLois Curfman McInnes b->rowners[0] = 0; 148944cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 149044cd7ae7SLois Curfman McInnes b->rowners[i] += b->rowners[i-1]; 14918a729477SBarry Smith } 149244cd7ae7SLois Curfman McInnes b->rstart = b->rowners[b->rank]; 149344cd7ae7SLois Curfman McInnes b->rend = b->rowners[b->rank+1]; 149444cd7ae7SLois Curfman McInnes MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm); 149544cd7ae7SLois Curfman McInnes b->cowners[0] = 0; 149644cd7ae7SLois Curfman McInnes for ( i=2; i<=b->size; i++ ) { 149744cd7ae7SLois Curfman McInnes b->cowners[i] += b->cowners[i-1]; 14988a729477SBarry Smith } 149944cd7ae7SLois Curfman McInnes b->cstart = b->cowners[b->rank]; 150044cd7ae7SLois Curfman McInnes b->cend = b->cowners[b->rank+1]; 15018a729477SBarry Smith 15025ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 150344cd7ae7SLois Curfman McInnes ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr); 150444cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->A); 15057b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 150644cd7ae7SLois Curfman McInnes ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr); 150744cd7ae7SLois Curfman McInnes PLogObjectParent(B,b->B); 15088a729477SBarry Smith 15091eb62cbbSBarry Smith /* build cache for off array entries formed */ 151044cd7ae7SLois Curfman McInnes ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr); 151144cd7ae7SLois Curfman McInnes b->colmap = 0; 151244cd7ae7SLois Curfman McInnes b->garray = 0; 151344cd7ae7SLois Curfman McInnes b->roworiented = 1; 15148a729477SBarry Smith 15151eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 151644cd7ae7SLois Curfman McInnes b->lvec = 0; 151744cd7ae7SLois Curfman McInnes b->Mvctx = 0; 15188a729477SBarry Smith 15197a0afa10SBarry Smith /* stuff for MatGetRow() */ 152044cd7ae7SLois Curfman McInnes b->rowindices = 0; 152144cd7ae7SLois Curfman McInnes b->rowvalues = 0; 152244cd7ae7SLois Curfman McInnes b->getrowactive = PETSC_FALSE; 15237a0afa10SBarry Smith 152444cd7ae7SLois Curfman McInnes *A = B; 1525d6dfbf8fSBarry Smith return 0; 1526d6dfbf8fSBarry Smith } 1527c74985f6SBarry Smith 15283d1612f7SBarry Smith static int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues) 1529d6dfbf8fSBarry Smith { 1530d6dfbf8fSBarry Smith Mat mat; 1531416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 1532a1b97e82SLois Curfman McInnes int ierr, len=0, flg; 1533d6dfbf8fSBarry Smith 1534416022c9SBarry Smith *newmat = 0; 15350452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm); 1536a5a9c739SBarry Smith PLogObjectCreate(mat); 15370452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1538416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 153944a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 154044a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1541d6dfbf8fSBarry Smith mat->factor = matin->factor; 1542c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1543d6dfbf8fSBarry Smith 154444cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 154544cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 154644cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 154744cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 1548d6dfbf8fSBarry Smith 1549416022c9SBarry Smith a->rstart = oldmat->rstart; 1550416022c9SBarry Smith a->rend = oldmat->rend; 1551416022c9SBarry Smith a->cstart = oldmat->cstart; 1552416022c9SBarry Smith a->cend = oldmat->cend; 155317699dbbSLois Curfman McInnes a->size = oldmat->size; 155417699dbbSLois Curfman McInnes a->rank = oldmat->rank; 155564119d58SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 1556bcd2baecSBarry Smith a->rowvalues = 0; 1557bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1558d6dfbf8fSBarry Smith 15590452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 156017699dbbSLois Curfman McInnes PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ)); 156117699dbbSLois Curfman McInnes PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 1562416022c9SBarry Smith ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 15632ee70a88SLois Curfman McInnes if (oldmat->colmap) { 15640452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1565416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1566416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1567416022c9SBarry Smith } else a->colmap = 0; 1568ec8511deSBarry Smith if (oldmat->garray && (len = ((Mat_SeqAIJ *) (oldmat->B->data))->n)) { 15690452661fSBarry Smith a->garray = (int *) PetscMalloc(len*sizeof(int)); CHKPTRQ(a->garray); 1570464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 1571416022c9SBarry Smith PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1572416022c9SBarry Smith } else a->garray = 0; 1573d6dfbf8fSBarry Smith 1574416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1575416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1576a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1577416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 1578416022c9SBarry Smith ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 1579416022c9SBarry Smith PLogObjectParent(mat,a->A); 1580416022c9SBarry Smith ierr = MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr); 1581416022c9SBarry Smith PLogObjectParent(mat,a->B); 15825dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 158325cdf11fSBarry Smith if (flg) { 1584682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1585682d7d0cSBarry Smith } 15868a729477SBarry Smith *newmat = mat; 15878a729477SBarry Smith return 0; 15888a729477SBarry Smith } 1589416022c9SBarry Smith 159077c4ece6SBarry Smith #include "sys.h" 1591416022c9SBarry Smith 159219bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat) 1593416022c9SBarry Smith { 1594d65a2f8fSBarry Smith Mat A; 1595416022c9SBarry Smith int i, nz, ierr, j,rstart, rend, fd; 1596d65a2f8fSBarry Smith Scalar *vals,*svals; 159719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 1598416022c9SBarry Smith MPI_Status status; 159917699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1600d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 160119bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 1602416022c9SBarry Smith 160317699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 160417699dbbSLois Curfman McInnes if (!rank) { 160590ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 160677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr); 1607c456f294SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIAIJ:not matrix object"); 1608416022c9SBarry Smith } 1609416022c9SBarry Smith 1610416022c9SBarry Smith MPI_Bcast(header+1,3,MPI_INT,0,comm); 1611416022c9SBarry Smith M = header[1]; N = header[2]; 1612416022c9SBarry Smith /* determine ownership of all rows */ 161317699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 16140452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1615416022c9SBarry Smith MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 1616416022c9SBarry Smith rowners[0] = 0; 161717699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1618416022c9SBarry Smith rowners[i] += rowners[i-1]; 1619416022c9SBarry Smith } 162017699dbbSLois Curfman McInnes rstart = rowners[rank]; 162117699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1622416022c9SBarry Smith 1623416022c9SBarry Smith /* distribute row lengths to all processors */ 16240452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1625416022c9SBarry Smith offlens = ourlens + (rend-rstart); 162617699dbbSLois Curfman McInnes if (!rank) { 16270452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 162877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 16290452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 163017699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1631d65a2f8fSBarry Smith MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 16320452661fSBarry Smith PetscFree(sndcounts); 1633416022c9SBarry Smith } 1634416022c9SBarry Smith else { 1635416022c9SBarry Smith MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 1636416022c9SBarry Smith } 1637416022c9SBarry Smith 163817699dbbSLois Curfman McInnes if (!rank) { 1639416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 16400452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1641cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 164217699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1643416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1644416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1645416022c9SBarry Smith } 1646416022c9SBarry Smith } 16470452661fSBarry Smith PetscFree(rowlengths); 1648416022c9SBarry Smith 1649416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1650416022c9SBarry Smith maxnz = 0; 165117699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 16520452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1653416022c9SBarry Smith } 16540452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1655416022c9SBarry Smith 1656416022c9SBarry Smith /* read in my part of the matrix column indices */ 1657416022c9SBarry Smith nz = procsnz[0]; 16580452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 165977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr); 1660d65a2f8fSBarry Smith 1661d65a2f8fSBarry Smith /* read in every one elses and ship off */ 166217699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1663d65a2f8fSBarry Smith nz = procsnz[i]; 166477c4ece6SBarry Smith ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 1665d65a2f8fSBarry Smith MPI_Send(cols,nz,MPI_INT,i,tag,comm); 1666d65a2f8fSBarry Smith } 16670452661fSBarry Smith PetscFree(cols); 1668416022c9SBarry Smith } 1669416022c9SBarry Smith else { 1670416022c9SBarry Smith /* determine buffer space needed for message */ 1671416022c9SBarry Smith nz = 0; 1672416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1673416022c9SBarry Smith nz += ourlens[i]; 1674416022c9SBarry Smith } 16750452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1676416022c9SBarry Smith 1677416022c9SBarry Smith /* receive message of column indices*/ 1678d65a2f8fSBarry Smith MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 1679416022c9SBarry Smith MPI_Get_count(&status,MPI_INT,&maxnz); 1680c456f294SBarry Smith if (maxnz != nz) SETERRQ(1,"MatLoad_MPIAIJ:something is wrong with file"); 1681416022c9SBarry Smith } 1682416022c9SBarry Smith 1683416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 1684cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 1685416022c9SBarry Smith jj = 0; 1686416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1687416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 1688d65a2f8fSBarry Smith if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1689416022c9SBarry Smith jj++; 1690416022c9SBarry Smith } 1691416022c9SBarry Smith } 1692d65a2f8fSBarry Smith 1693d65a2f8fSBarry Smith /* create our matrix */ 1694416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1695416022c9SBarry Smith ourlens[i] -= offlens[i]; 1696416022c9SBarry Smith } 1697d65a2f8fSBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 1698d65a2f8fSBarry Smith A = *newmat; 1699d65a2f8fSBarry Smith MatSetOption(A,COLUMNS_SORTED); 1700d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1701d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 1702d65a2f8fSBarry Smith } 1703416022c9SBarry Smith 170417699dbbSLois Curfman McInnes if (!rank) { 17050452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1706416022c9SBarry Smith 1707416022c9SBarry Smith /* read in my part of the matrix numerical values */ 1708416022c9SBarry Smith nz = procsnz[0]; 170977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 1710d65a2f8fSBarry Smith 1711d65a2f8fSBarry Smith /* insert into matrix */ 1712d65a2f8fSBarry Smith jj = rstart; 1713d65a2f8fSBarry Smith smycols = mycols; 1714d65a2f8fSBarry Smith svals = vals; 1715d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1716d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1717d65a2f8fSBarry Smith smycols += ourlens[i]; 1718d65a2f8fSBarry Smith svals += ourlens[i]; 1719d65a2f8fSBarry Smith jj++; 1720416022c9SBarry Smith } 1721416022c9SBarry Smith 1722d65a2f8fSBarry Smith /* read in other processors and ship out */ 172317699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1724416022c9SBarry Smith nz = procsnz[i]; 172577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 1726d65a2f8fSBarry Smith MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 1727416022c9SBarry Smith } 17280452661fSBarry Smith PetscFree(procsnz); 1729416022c9SBarry Smith } 1730d65a2f8fSBarry Smith else { 1731d65a2f8fSBarry Smith /* receive numeric values */ 17320452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1733416022c9SBarry Smith 1734d65a2f8fSBarry Smith /* receive message of values*/ 1735d65a2f8fSBarry Smith MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 1736d65a2f8fSBarry Smith MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 1737c456f294SBarry Smith if (maxnz != nz) SETERRQ(1,"MatLoad_MPIAIJ:something is wrong with file"); 1738d65a2f8fSBarry Smith 1739d65a2f8fSBarry Smith /* insert into matrix */ 1740d65a2f8fSBarry Smith jj = rstart; 1741d65a2f8fSBarry Smith smycols = mycols; 1742d65a2f8fSBarry Smith svals = vals; 1743d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1744d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1745d65a2f8fSBarry Smith smycols += ourlens[i]; 1746d65a2f8fSBarry Smith svals += ourlens[i]; 1747d65a2f8fSBarry Smith jj++; 1748d65a2f8fSBarry Smith } 1749d65a2f8fSBarry Smith } 17500452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1751d65a2f8fSBarry Smith 1752d65a2f8fSBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 1753d65a2f8fSBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 1754416022c9SBarry Smith return 0; 1755416022c9SBarry Smith } 1756