1cb512458SBarry Smith #ifndef lint 2*bcd2baecSBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.130 1996/03/04 05:15:55 bsmith Exp bsmith $"; 3cb512458SBarry Smith #endif 48a729477SBarry Smith 51eb62cbbSBarry Smith #include "mpiaij.h" 68a729477SBarry Smith #include "vec/vecimpl.h" 7d6dfbf8fSBarry Smith #include "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 32778b31e54SBarry Smith ierr = ISGetLocalSize(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); 443416022c9SBarry Smith ierr = MatMult(a->A,xx,yy); CHKERRQ(ierr); 44464119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr); 445416022c9SBarry Smith ierr = MatMultAdd(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); 454416022c9SBarry Smith ierr = MatMultAdd(a->A,xx,yy,zz); CHKERRQ(ierr); 45564119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr); 456416022c9SBarry Smith ierr = MatMultAdd(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 */ 466416022c9SBarry Smith ierr = MatMultTrans(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 */ 471416022c9SBarry Smith ierr = MatMultTrans(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 */ 486416022c9SBarry Smith ierr = MatMultTrans(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 */ 491416022c9SBarry Smith ierr = MatMultTransAdd(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; 507416022c9SBarry Smith return MatGetDiagonal(a->A,v); 5081eb62cbbSBarry Smith } 5091eb62cbbSBarry Smith 510052efed2SBarry Smith static int MatScale_MPIAIJ(Scalar *aa,Mat A) 511052efed2SBarry Smith { 512052efed2SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 513052efed2SBarry Smith int ierr; 514052efed2SBarry Smith ierr = MatScale(aa,a->A); CHKERRQ(ierr); 515052efed2SBarry Smith ierr = MatScale(aa,a->B); CHKERRQ(ierr); 516052efed2SBarry Smith return 0; 517052efed2SBarry Smith } 518052efed2SBarry Smith 51944a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj) 5201eb62cbbSBarry Smith { 5211eb62cbbSBarry Smith Mat mat = (Mat) obj; 52244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 5231eb62cbbSBarry Smith int ierr; 524a5a9c739SBarry Smith #if defined(PETSC_LOG) 5256e35fa57SLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N); 526a5a9c739SBarry Smith #endif 5270452661fSBarry Smith PetscFree(aij->rowners); 52878b31e54SBarry Smith ierr = MatDestroy(aij->A); CHKERRQ(ierr); 52978b31e54SBarry Smith ierr = MatDestroy(aij->B); CHKERRQ(ierr); 5300452661fSBarry Smith if (aij->colmap) PetscFree(aij->colmap); 5310452661fSBarry Smith if (aij->garray) PetscFree(aij->garray); 5321eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 533a56f8943SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 5347a0afa10SBarry Smith if (aij->rowvalues) PetscFree(aij->rowvalues); 5350452661fSBarry Smith PetscFree(aij); 536a5a9c739SBarry Smith PLogObjectDestroy(mat); 5370452661fSBarry Smith PetscHeaderDestroy(mat); 5381eb62cbbSBarry Smith return 0; 5391eb62cbbSBarry Smith } 5407857610eSBarry Smith #include "draw.h" 541b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 542ee50ffe9SBarry Smith 543416022c9SBarry Smith static int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer) 5441eb62cbbSBarry Smith { 545416022c9SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 546416022c9SBarry Smith int ierr; 547416022c9SBarry Smith 54817699dbbSLois Curfman McInnes if (aij->size == 1) { 549416022c9SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 550416022c9SBarry Smith } 551a523beb4SLois Curfman McInnes else SETERRQ(1,"MatView_MPIAIJ_Binary:Only uniprocessor output supported"); 552416022c9SBarry Smith return 0; 553416022c9SBarry Smith } 554416022c9SBarry Smith 555d7e8b826SBarry Smith static int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer) 556416022c9SBarry Smith { 55744a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 558dbb450caSBarry Smith Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data; 559a56f8943SBarry Smith int ierr, format,shift = C->indexshift,rank; 560d13ab20cSBarry Smith PetscObject vobj = (PetscObject) viewer; 561d636dbe3SBarry Smith FILE *fd; 562416022c9SBarry Smith 563416022c9SBarry Smith if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 564416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 56508480c60SBarry Smith if (format == FILE_FORMAT_INFO_DETAILED) { 56695e01e2fSLois Curfman McInnes int nz, nzalloc, mem, flg; 567a56f8943SBarry Smith MPI_Comm_rank(mat->comm,&rank); 568227d817aSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 569a56f8943SBarry Smith ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); 57095e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 571a56f8943SBarry Smith MPIU_Seq_begin(mat->comm,1); 57295e01e2fSLois Curfman McInnes if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n", 57395e01e2fSLois Curfman McInnes rank,aij->m,nz,nzalloc,mem); 57495e01e2fSLois Curfman McInnes else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n", 57595e01e2fSLois Curfman McInnes rank,aij->m,nz,nzalloc,mem); 57608480c60SBarry Smith ierr = MatGetInfo(aij->A,MAT_LOCAL,&nz,&nzalloc,&mem); 57795e01e2fSLois Curfman McInnes fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,nz); 57808480c60SBarry Smith ierr = MatGetInfo(aij->B,MAT_LOCAL,&nz,&nzalloc,&mem); 57995e01e2fSLois Curfman McInnes fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,nz); 580a56f8943SBarry Smith fflush(fd); 581a56f8943SBarry Smith MPIU_Seq_end(mat->comm,1); 582a40aa06bSLois Curfman McInnes ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr); 583416022c9SBarry Smith return 0; 584416022c9SBarry Smith } 58508480c60SBarry Smith else if (format == FILE_FORMAT_INFO) { 58608480c60SBarry Smith return 0; 58708480c60SBarry Smith } 588416022c9SBarry Smith } 589416022c9SBarry Smith 590416022c9SBarry Smith if (vobj->type == ASCII_FILE_VIEWER) { 591227d817aSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 5927f813858SBarry Smith MPIU_Seq_begin(mat->comm,1); 593d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 59417699dbbSLois Curfman McInnes aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 5951eb62cbbSBarry Smith aij->cend); 59678b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 59778b31e54SBarry Smith ierr = MatView(aij->B,viewer); CHKERRQ(ierr); 598d13ab20cSBarry Smith fflush(fd); 5997f813858SBarry Smith MPIU_Seq_end(mat->comm,1); 600d13ab20cSBarry Smith } 601416022c9SBarry Smith else { 602a56f8943SBarry Smith int size = aij->size; 603a56f8943SBarry Smith rank = aij->rank; 60417699dbbSLois Curfman McInnes if (size == 1) { 60578b31e54SBarry Smith ierr = MatView(aij->A,viewer); CHKERRQ(ierr); 60695373324SBarry Smith } 60795373324SBarry Smith else { 60895373324SBarry Smith /* assemble the entire matrix onto first processor. */ 60995373324SBarry Smith Mat A; 610ec8511deSBarry Smith Mat_SeqAIJ *Aloc; 6112eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 61295373324SBarry Smith Scalar *a; 6132ee70a88SLois Curfman McInnes 61417699dbbSLois Curfman McInnes if (!rank) { 615b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 616c451ab25SLois Curfman McInnes CHKERRQ(ierr); 61795373324SBarry Smith } 61895373324SBarry Smith else { 619b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A); 620c451ab25SLois Curfman McInnes CHKERRQ(ierr); 62195373324SBarry Smith } 622464493b3SBarry Smith PLogObjectParent(mat,A); 623416022c9SBarry Smith 62495373324SBarry Smith /* copy over the A part */ 625ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->A->data; 6262ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 62795373324SBarry Smith row = aij->rstart; 628dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;} 62995373324SBarry Smith for ( i=0; i<m; i++ ) { 630416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr); 63195373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 63295373324SBarry Smith } 6332ee70a88SLois Curfman McInnes aj = Aloc->j; 634dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;} 63595373324SBarry Smith 63695373324SBarry Smith /* copy over the B part */ 637ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) aij->B->data; 6382ee70a88SLois Curfman McInnes m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a; 63995373324SBarry Smith row = aij->rstart; 6400452661fSBarry Smith ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols); 641dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];} 64295373324SBarry Smith for ( i=0; i<m; i++ ) { 643416022c9SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr); 64495373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 64595373324SBarry Smith } 6460452661fSBarry Smith PetscFree(ct); 64778b31e54SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 64878b31e54SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 64917699dbbSLois Curfman McInnes if (!rank) { 65078b31e54SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr); 65195373324SBarry Smith } 65278b31e54SBarry Smith ierr = MatDestroy(A); CHKERRQ(ierr); 65395373324SBarry Smith } 65495373324SBarry Smith } 6551eb62cbbSBarry Smith return 0; 6561eb62cbbSBarry Smith } 6571eb62cbbSBarry Smith 658416022c9SBarry Smith static int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 659416022c9SBarry Smith { 660416022c9SBarry Smith Mat mat = (Mat) obj; 661416022c9SBarry Smith int ierr; 662416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 663416022c9SBarry Smith 664416022c9SBarry Smith if (!viewer) { 665416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 666416022c9SBarry Smith } 667416022c9SBarry Smith if (vobj->cookie == DRAW_COOKIE && vobj->type == NULLWINDOW) return 0; 668416022c9SBarry Smith else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILE_VIEWER) { 669d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 670416022c9SBarry Smith } 671416022c9SBarry Smith else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILES_VIEWER) { 672d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 673d7e8b826SBarry Smith } 674d7e8b826SBarry Smith else if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 675d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 676416022c9SBarry Smith } 677416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 678d7e8b826SBarry Smith ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr); 679416022c9SBarry Smith } 680416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 681416022c9SBarry Smith return MatView_MPIAIJ_Binary(mat,viewer); 682416022c9SBarry Smith } 683416022c9SBarry Smith return 0; 684416022c9SBarry Smith } 685416022c9SBarry Smith 686ec8511deSBarry Smith extern int MatMarkDiag_SeqAIJ(Mat); 6871eb62cbbSBarry Smith /* 6881eb62cbbSBarry Smith This has to provide several versions. 6891eb62cbbSBarry Smith 6901eb62cbbSBarry Smith 1) per sequential 6911eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 6921eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 693d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 6941eb62cbbSBarry Smith */ 695fc76ce87SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag, 696dbb450caSBarry Smith double fshift,int its,Vec xx) 6978a729477SBarry Smith { 69844a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 699d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 700ec8511deSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data; 7016abc6512SBarry Smith Scalar zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts; 7026abc6512SBarry Smith int ierr,*idx, *diag; 703416022c9SBarry Smith int n = mat->n, m = mat->m, i,shift = A->indexshift; 704da3a660dSBarry Smith Vec tt; 7058a729477SBarry Smith 706d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 707dbb450caSBarry Smith xs = x + shift; /* shift by one for index start of 1 */ 708dbb450caSBarry Smith ls = ls + shift; 709ec8511deSBarry Smith if (!A->diag) {if ((ierr = MatMarkDiag_SeqAIJ(AA))) return ierr;} 710d6dfbf8fSBarry Smith diag = A->diag; 711acb40c82SBarry Smith if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) { 71248d91487SBarry Smith SETERRQ(1,"MatRelax_MPIAIJ:Option not supported"); 713acb40c82SBarry Smith } 714da3a660dSBarry Smith if (flag & SOR_EISENSTAT) { 715da3a660dSBarry Smith /* Let A = L + U + D; where L is lower trianglar, 716da3a660dSBarry Smith U is upper triangular, E is diagonal; This routine applies 717da3a660dSBarry Smith 718da3a660dSBarry Smith (L + E)^{-1} A (U + E)^{-1} 719da3a660dSBarry Smith 720da3a660dSBarry Smith to a vector efficiently using Eisenstat's trick. This is for 721da3a660dSBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 722da3a660dSBarry Smith is the relaxation factor. 723da3a660dSBarry Smith */ 72478b31e54SBarry Smith ierr = VecDuplicate(xx,&tt); CHKERRQ(ierr); 725464493b3SBarry Smith PLogObjectParent(matin,tt); 726da3a660dSBarry Smith VecGetArray(tt,&t); 727da3a660dSBarry Smith scale = (2.0/omega) - 1.0; 728da3a660dSBarry Smith /* x = (E + U)^{-1} b */ 729da3a660dSBarry Smith VecSet(&zero,mat->lvec); 73064119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 73178b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 732da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 733da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 734dbb450caSBarry Smith idx = A->j + diag[i] + !shift; 735dbb450caSBarry Smith v = A->a + diag[i] + !shift; 736da3a660dSBarry Smith sum = b[i]; 737da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 738dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 739da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 740dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 741dbb450caSBarry Smith v = B->a + B->i[i] + shift; 742da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 743da3a660dSBarry Smith x[i] = omega*(sum/d); 744da3a660dSBarry Smith } 74564119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 74678b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 747da3a660dSBarry Smith 748da3a660dSBarry Smith /* t = b - (2*E - D)x */ 749da3a660dSBarry Smith v = A->a; 750dbb450caSBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 751da3a660dSBarry Smith 752da3a660dSBarry Smith /* t = (E + L)^{-1}t */ 753dbb450caSBarry Smith ts = t + shift; /* shifted by one for index start of a or mat->j*/ 754da3a660dSBarry Smith diag = A->diag; 755da3a660dSBarry Smith VecSet(&zero,mat->lvec); 75664119d58SLois Curfman McInnes ierr = VecPipelineBegin(tt,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 75778b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 758da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 759da3a660dSBarry Smith n = diag[i] - A->i[i]; 760dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 761dbb450caSBarry Smith v = A->a + A->i[i] + shift; 762da3a660dSBarry Smith sum = t[i]; 763da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 764dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 765da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 766dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 767dbb450caSBarry Smith v = B->a + B->i[i] + shift; 768da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 769da3a660dSBarry Smith t[i] = omega*(sum/d); 770da3a660dSBarry Smith } 77164119d58SLois Curfman McInnes ierr = VecPipelineEnd(tt,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 77278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 773da3a660dSBarry Smith /* x = x + t */ 774da3a660dSBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 775da3a660dSBarry Smith VecDestroy(tt); 776da3a660dSBarry Smith return 0; 777da3a660dSBarry Smith } 778da3a660dSBarry Smith 7791eb62cbbSBarry Smith 780d6dfbf8fSBarry Smith if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){ 781da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 782da3a660dSBarry Smith VecSet(&zero,mat->lvec); VecSet(&zero,xx); 783da3a660dSBarry Smith } 784da3a660dSBarry Smith else { 78564119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 78678b31e54SBarry Smith CHKERRQ(ierr); 78764119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 78878b31e54SBarry Smith CHKERRQ(ierr); 789da3a660dSBarry Smith } 790d6dfbf8fSBarry Smith while (its--) { 791d6dfbf8fSBarry Smith /* go down through the rows */ 79264119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 79378b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 794d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 795d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 796dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 797dbb450caSBarry Smith v = A->a + A->i[i] + shift; 798d6dfbf8fSBarry Smith sum = b[i]; 799d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 800dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 801d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 802dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 803dbb450caSBarry Smith v = B->a + B->i[i] + shift; 804d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 805d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 806d6dfbf8fSBarry Smith } 80764119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 80878b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 809d6dfbf8fSBarry Smith /* come up through the rows */ 81064119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 81178b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 812d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 813d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 814dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 815dbb450caSBarry Smith v = A->a + A->i[i] + shift; 816d6dfbf8fSBarry Smith sum = b[i]; 817d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 818dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 819d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 820dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 821dbb450caSBarry Smith v = B->a + B->i[i] + shift; 822d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 823d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 824d6dfbf8fSBarry Smith } 82564119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 82678b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 827d6dfbf8fSBarry Smith } 828d6dfbf8fSBarry Smith } 829d6dfbf8fSBarry Smith else if (flag & SOR_FORWARD_SWEEP){ 830da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 831da3a660dSBarry Smith VecSet(&zero,mat->lvec); 83264119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 83378b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 834da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 835da3a660dSBarry Smith n = diag[i] - A->i[i]; 836dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 837dbb450caSBarry Smith v = A->a + A->i[i] + shift; 838da3a660dSBarry Smith sum = b[i]; 839da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 840dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 841da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 842dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 843dbb450caSBarry Smith v = B->a + B->i[i] + shift; 844da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 845da3a660dSBarry Smith x[i] = omega*(sum/d); 846da3a660dSBarry Smith } 84764119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 84878b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 849da3a660dSBarry Smith its--; 850da3a660dSBarry Smith } 851d6dfbf8fSBarry Smith while (its--) { 85264119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 85378b31e54SBarry Smith CHKERRQ(ierr); 85464119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx); 85578b31e54SBarry Smith CHKERRQ(ierr); 85664119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 85778b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 858d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 859d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 860dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 861dbb450caSBarry Smith v = A->a + A->i[i] + shift; 862d6dfbf8fSBarry Smith sum = b[i]; 863d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 864dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 865d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 866dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 867dbb450caSBarry Smith v = B->a + B->i[i] + shift; 868d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 869d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 870d6dfbf8fSBarry Smith } 87164119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN, 87278b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 873d6dfbf8fSBarry Smith } 874d6dfbf8fSBarry Smith } 875d6dfbf8fSBarry Smith else if (flag & SOR_BACKWARD_SWEEP){ 876da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 877da3a660dSBarry Smith VecSet(&zero,mat->lvec); 87864119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 87978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 880da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 881da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 882dbb450caSBarry Smith idx = A->j + diag[i] + !shift; 883dbb450caSBarry Smith v = A->a + diag[i] + !shift; 884da3a660dSBarry Smith sum = b[i]; 885da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 886dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 887da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 888dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 889dbb450caSBarry Smith v = B->a + B->i[i] + shift; 890da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 891da3a660dSBarry Smith x[i] = omega*(sum/d); 892da3a660dSBarry Smith } 89364119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 89478b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 895da3a660dSBarry Smith its--; 896da3a660dSBarry Smith } 897d6dfbf8fSBarry Smith while (its--) { 89864119d58SLois Curfman McInnes ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_DOWN, 89978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 90064119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_DOWN, 90178b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 90264119d58SLois Curfman McInnes ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 90378b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 904d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 905d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 906dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 907dbb450caSBarry Smith v = A->a + A->i[i] + shift; 908d6dfbf8fSBarry Smith sum = b[i]; 909d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 910dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 911d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 912dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 913dbb450caSBarry Smith v = B->a + B->i[i] + shift; 914d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 915d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 916d6dfbf8fSBarry Smith } 91764119d58SLois Curfman McInnes ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP, 91878b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 919d6dfbf8fSBarry Smith } 920d6dfbf8fSBarry Smith } 921d6dfbf8fSBarry Smith else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 922da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 923dbb450caSBarry Smith return MatRelax(mat->A,bb,omega,flag,fshift,its,xx); 924da3a660dSBarry Smith } 92564119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 92678b31e54SBarry Smith CHKERRQ(ierr); 92764119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 92878b31e54SBarry Smith CHKERRQ(ierr); 929d6dfbf8fSBarry Smith while (its--) { 930d6dfbf8fSBarry Smith /* go down through the rows */ 931d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 932d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 933dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 934dbb450caSBarry Smith v = A->a + A->i[i] + shift; 935d6dfbf8fSBarry Smith sum = b[i]; 936d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 937dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 938d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 939dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 940dbb450caSBarry Smith v = B->a + B->i[i] + shift; 941d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 942d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 943d6dfbf8fSBarry Smith } 944d6dfbf8fSBarry Smith /* come up through the rows */ 945d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 946d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 947dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 948dbb450caSBarry Smith v = A->a + A->i[i] + shift; 949d6dfbf8fSBarry Smith sum = b[i]; 950d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 951dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 952d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 953dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 954dbb450caSBarry Smith v = B->a + B->i[i] + shift; 955d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 956d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 957d6dfbf8fSBarry Smith } 958d6dfbf8fSBarry Smith } 959d6dfbf8fSBarry Smith } 960d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 961da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 962dbb450caSBarry Smith return MatRelax(mat->A,bb,omega,flag,fshift,its,xx); 963da3a660dSBarry Smith } 96464119d58SLois Curfman McInnes ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 96578b31e54SBarry Smith CHKERRQ(ierr); 96664119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx); 96778b31e54SBarry Smith CHKERRQ(ierr); 968d6dfbf8fSBarry Smith while (its--) { 969d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 970d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 971dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 972dbb450caSBarry Smith v = A->a + A->i[i] + shift; 973d6dfbf8fSBarry Smith sum = b[i]; 974d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 975dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 976d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 977dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 978dbb450caSBarry Smith v = B->a + B->i[i] + shift; 979d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 980d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 981d6dfbf8fSBarry Smith } 982d6dfbf8fSBarry Smith } 983d6dfbf8fSBarry Smith } 984d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 985da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 986dbb450caSBarry Smith return MatRelax(mat->A,bb,omega,flag,fshift,its,xx); 987da3a660dSBarry Smith } 98864119d58SLois Curfman McInnes ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL, 98978b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 99064119d58SLois Curfman McInnes ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL, 99178b31e54SBarry Smith mat->Mvctx); CHKERRQ(ierr); 992d6dfbf8fSBarry Smith while (its--) { 993d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 994d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 995dbb450caSBarry Smith idx = A->j + A->i[i] + shift; 996dbb450caSBarry Smith v = A->a + A->i[i] + shift; 997d6dfbf8fSBarry Smith sum = b[i]; 998d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 999dbb450caSBarry Smith d = fshift + A->a[diag[i]+shift]; 1000d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 1001dbb450caSBarry Smith idx = B->j + B->i[i] + shift; 1002dbb450caSBarry Smith v = B->a + B->i[i] + shift; 1003d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 1004d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 1005d6dfbf8fSBarry Smith } 1006d6dfbf8fSBarry Smith } 1007d6dfbf8fSBarry Smith } 10088a729477SBarry Smith return 0; 10098a729477SBarry Smith } 1010a66be287SLois Curfman McInnes 1011fc76ce87SLois Curfman McInnes static int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,int *nz, 1012fc76ce87SLois Curfman McInnes int *nzalloc,int *mem) 1013a66be287SLois Curfman McInnes { 1014a66be287SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1015a66be287SLois Curfman McInnes Mat A = mat->A, B = mat->B; 1016a66be287SLois Curfman McInnes int ierr, isend[3], irecv[3], nzA, nzallocA, memA; 1017a66be287SLois Curfman McInnes 101878b31e54SBarry Smith ierr = MatGetInfo(A,MAT_LOCAL,&nzA,&nzallocA,&memA); CHKERRQ(ierr); 101978b31e54SBarry Smith ierr = MatGetInfo(B,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr); 1020a66be287SLois Curfman McInnes isend[0] += nzA; isend[1] += nzallocA; isend[2] += memA; 1021a66be287SLois Curfman McInnes if (flag == MAT_LOCAL) { 1022*bcd2baecSBarry Smith if (nz) *nz = isend[0]; 1023*bcd2baecSBarry Smith if (nzalloc) *nzalloc = isend[1]; 1024*bcd2baecSBarry Smith if (mem) *mem = isend[2]; 1025a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 1026d65a2f8fSBarry Smith MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_MAX,matin->comm); 1027*bcd2baecSBarry Smith if (nz) *nz = irecv[0]; 1028*bcd2baecSBarry Smith if (nzalloc) *nzalloc = irecv[1]; 1029*bcd2baecSBarry Smith if (mem) *mem = irecv[2]; 1030a66be287SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 1031d65a2f8fSBarry Smith MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_SUM,matin->comm); 1032*bcd2baecSBarry Smith if (nz) *nz = irecv[0]; 1033*bcd2baecSBarry Smith if (nzalloc) *nzalloc = irecv[1]; 1034*bcd2baecSBarry Smith if (mem) *mem = irecv[2]; 1035a66be287SLois Curfman McInnes } 1036a66be287SLois Curfman McInnes return 0; 1037a66be287SLois Curfman McInnes } 1038a66be287SLois Curfman McInnes 1039299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*); 1040299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*); 1041299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double); 1042299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *); 1043299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec); 1044299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec); 1045299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec); 1046299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec); 1047299609e3SLois Curfman McInnes 1048416022c9SBarry Smith static int MatSetOption_MPIAIJ(Mat A,MatOption op) 1049c74985f6SBarry Smith { 1050c0bbcb79SLois Curfman McInnes Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1051c74985f6SBarry Smith 1052c0bbcb79SLois Curfman McInnes if (op == NO_NEW_NONZERO_LOCATIONS || 1053c0bbcb79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 1054c0bbcb79SLois Curfman McInnes op == COLUMNS_SORTED || 1055c0bbcb79SLois Curfman McInnes op == ROW_ORIENTED) { 1056c0bbcb79SLois Curfman McInnes MatSetOption(a->A,op); 1057c0bbcb79SLois Curfman McInnes MatSetOption(a->B,op); 1058c74985f6SBarry Smith } 1059c0bbcb79SLois Curfman McInnes else if (op == ROWS_SORTED || 1060c0bbcb79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 1061c0bbcb79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 1062c0bbcb79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 1063c0bbcb79SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_MPIAIJ:Option ignored\n"); 10644b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) { 10654b0e389bSBarry Smith a->roworiented = 0; 10664b0e389bSBarry Smith MatSetOption(a->A,op); 10674b0e389bSBarry Smith MatSetOption(a->B,op); 10684b0e389bSBarry Smith } 1069c0bbcb79SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 10704d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:NO_NEW_DIAGONALS");} 1071c0bbcb79SLois Curfman McInnes else 10724d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:unknown option");} 1073c74985f6SBarry Smith return 0; 1074c74985f6SBarry Smith } 1075c74985f6SBarry Smith 1076d1710a03SLois Curfman McInnes static int MatGetSize_MPIAIJ(Mat matin,int *m,int *n) 1077c74985f6SBarry Smith { 107844a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1079c74985f6SBarry Smith *m = mat->M; *n = mat->N; 1080c74985f6SBarry Smith return 0; 1081c74985f6SBarry Smith } 1082c74985f6SBarry Smith 1083d1710a03SLois Curfman McInnes static int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n) 1084c74985f6SBarry Smith { 108544a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1086b7c46309SBarry Smith *m = mat->m; *n = mat->N; 1087c74985f6SBarry Smith return 0; 1088c74985f6SBarry Smith } 1089c74985f6SBarry Smith 1090d1710a03SLois Curfman McInnes static int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n) 1091c74985f6SBarry Smith { 109244a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 1093c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 1094c74985f6SBarry Smith return 0; 1095c74985f6SBarry Smith } 1096c74985f6SBarry Smith 10976d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 10986d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**); 10996d84be18SBarry Smith 11006d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 110139e00950SLois Curfman McInnes { 1102154123eaSLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 110370f0671dSBarry Smith Scalar *vworkA, *vworkB, **pvA, **pvB,*v_p; 1104154123eaSLois Curfman McInnes int i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart; 1105154123eaSLois Curfman McInnes int nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend; 110670f0671dSBarry Smith int *cmap, *idx_p; 110739e00950SLois Curfman McInnes 11087a0afa10SBarry Smith if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,"MatGetRow_MPIAIJ:Already active"); 11097a0afa10SBarry Smith mat->getrowactive = PETSC_TRUE; 11107a0afa10SBarry Smith 111170f0671dSBarry Smith if (!mat->rowvalues && (idx || v)) { 11127a0afa10SBarry Smith /* 11137a0afa10SBarry Smith allocate enough space to hold information from the longest row. 11147a0afa10SBarry Smith */ 11157a0afa10SBarry Smith Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data; 11167a0afa10SBarry Smith int max = 1,n = mat->n,tmp; 11177a0afa10SBarry Smith for ( i=0; i<n; i++ ) { 11187a0afa10SBarry Smith tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; 11197a0afa10SBarry Smith if (max < tmp) { max = tmp; } 11207a0afa10SBarry Smith } 11217a0afa10SBarry Smith mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar))); 11227a0afa10SBarry Smith CHKPTRQ(mat->rowvalues); 11237a0afa10SBarry Smith mat->rowindices = (int *) (mat->rowvalues + max); 11247a0afa10SBarry Smith } 11257a0afa10SBarry Smith 11267a0afa10SBarry Smith 1127b49de8d1SLois Curfman McInnes if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIAIJ:Only local rows") 1128abc0e9e4SLois Curfman McInnes lrow = row - rstart; 112939e00950SLois Curfman McInnes 1130154123eaSLois Curfman McInnes pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB; 1131154123eaSLois Curfman McInnes if (!v) {pvA = 0; pvB = 0;} 1132154123eaSLois Curfman McInnes if (!idx) {pcA = 0; if (!v) pcB = 0;} 11336d84be18SBarry Smith ierr = MatGetRow_SeqAIJ(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 11346d84be18SBarry Smith ierr = MatGetRow_SeqAIJ(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 1135154123eaSLois Curfman McInnes nztot = nzA + nzB; 1136154123eaSLois Curfman McInnes 113770f0671dSBarry Smith cmap = mat->garray; 1138154123eaSLois Curfman McInnes if (v || idx) { 1139154123eaSLois Curfman McInnes if (nztot) { 1140154123eaSLois Curfman McInnes /* Sort by increasing column numbers, assuming A and B already sorted */ 114170f0671dSBarry Smith int imark = -1; 1142154123eaSLois Curfman McInnes if (v) { 114370f0671dSBarry Smith *v = v_p = mat->rowvalues; 114439e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 114570f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) v_p[i] = vworkB[i]; 1146154123eaSLois Curfman McInnes else break; 1147154123eaSLois Curfman McInnes } 1148154123eaSLois Curfman McInnes imark = i; 114970f0671dSBarry Smith for ( i=0; i<nzA; i++ ) v_p[imark+i] = vworkA[i]; 115070f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) v_p[nzA+i] = vworkB[i]; 1151154123eaSLois Curfman McInnes } 1152154123eaSLois Curfman McInnes if (idx) { 115370f0671dSBarry Smith *idx = idx_p = mat->rowindices; 115470f0671dSBarry Smith if (imark > -1) { 115570f0671dSBarry Smith for ( i=0; i<imark; i++ ) { 115670f0671dSBarry Smith idx_p[i] = cmap[cworkB[i]]; 115770f0671dSBarry Smith } 115870f0671dSBarry Smith } else { 1159154123eaSLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 116070f0671dSBarry Smith if (cmap[cworkB[i]] < cstart) idx_p[i] = cmap[cworkB[i]]; 1161154123eaSLois Curfman McInnes else break; 1162154123eaSLois Curfman McInnes } 1163154123eaSLois Curfman McInnes imark = i; 116470f0671dSBarry Smith } 116570f0671dSBarry Smith for ( i=0; i<nzA; i++ ) idx_p[imark+i] = cstart + cworkA[i]; 116670f0671dSBarry Smith for ( i=imark; i<nzB; i++ ) idx_p[nzA+i] = cmap[cworkB[i]]; 116739e00950SLois Curfman McInnes } 116839e00950SLois Curfman McInnes } 116939e00950SLois Curfman McInnes else {*idx = 0; *v=0;} 1170154123eaSLois Curfman McInnes } 117139e00950SLois Curfman McInnes *nz = nztot; 11726d84be18SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr); 11736d84be18SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr); 117439e00950SLois Curfman McInnes return 0; 117539e00950SLois Curfman McInnes } 117639e00950SLois Curfman McInnes 11776d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 117839e00950SLois Curfman McInnes { 11797a0afa10SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 11807a0afa10SBarry Smith if (aij->getrowactive == PETSC_FALSE) { 11817a0afa10SBarry Smith SETERRQ(1,"MatRestoreRow_MPIAIJ:MatGetRow not called"); 11827a0afa10SBarry Smith } 11837a0afa10SBarry Smith aij->getrowactive = PETSC_FALSE; 118439e00950SLois Curfman McInnes return 0; 118539e00950SLois Curfman McInnes } 118639e00950SLois Curfman McInnes 1187cddf8d76SBarry Smith static int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm) 1188855ac2c5SLois Curfman McInnes { 1189855ac2c5SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 1190ec8511deSBarry Smith Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data; 1191416022c9SBarry Smith int ierr, i, j, cstart = aij->cstart,shift = amat->indexshift; 1192416022c9SBarry Smith double sum = 0.0; 119304ca555eSLois Curfman McInnes Scalar *v; 119404ca555eSLois Curfman McInnes 119517699dbbSLois Curfman McInnes if (aij->size == 1) { 119614183eadSLois Curfman McInnes ierr = MatNorm(aij->A,type,norm); CHKERRQ(ierr); 119737fa93a5SLois Curfman McInnes } else { 119804ca555eSLois Curfman McInnes if (type == NORM_FROBENIUS) { 119904ca555eSLois Curfman McInnes v = amat->a; 120004ca555eSLois Curfman McInnes for (i=0; i<amat->nz; i++ ) { 120104ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 120204ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 120304ca555eSLois Curfman McInnes #else 120404ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 120504ca555eSLois Curfman McInnes #endif 120604ca555eSLois Curfman McInnes } 120704ca555eSLois Curfman McInnes v = bmat->a; 120804ca555eSLois Curfman McInnes for (i=0; i<bmat->nz; i++ ) { 120904ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX) 121004ca555eSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 121104ca555eSLois Curfman McInnes #else 121204ca555eSLois Curfman McInnes sum += (*v)*(*v); v++; 121304ca555eSLois Curfman McInnes #endif 121404ca555eSLois Curfman McInnes } 12156d84be18SBarry Smith MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm); 121604ca555eSLois Curfman McInnes *norm = sqrt(*norm); 121704ca555eSLois Curfman McInnes } 121804ca555eSLois Curfman McInnes else if (type == NORM_1) { /* max column norm */ 121904ca555eSLois Curfman McInnes double *tmp, *tmp2; 122004ca555eSLois Curfman McInnes int *jj, *garray = aij->garray; 12210452661fSBarry Smith tmp = (double *) PetscMalloc( aij->N*sizeof(double) ); CHKPTRQ(tmp); 12220452661fSBarry Smith tmp2 = (double *) PetscMalloc( aij->N*sizeof(double) ); CHKPTRQ(tmp2); 1223cddf8d76SBarry Smith PetscMemzero(tmp,aij->N*sizeof(double)); 122404ca555eSLois Curfman McInnes *norm = 0.0; 122504ca555eSLois Curfman McInnes v = amat->a; jj = amat->j; 122604ca555eSLois Curfman McInnes for ( j=0; j<amat->nz; j++ ) { 1227579c6b6fSBarry Smith tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v); v++; 122804ca555eSLois Curfman McInnes } 122904ca555eSLois Curfman McInnes v = bmat->a; jj = bmat->j; 123004ca555eSLois Curfman McInnes for ( j=0; j<bmat->nz; j++ ) { 1231579c6b6fSBarry Smith tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++; 123204ca555eSLois Curfman McInnes } 12336d84be18SBarry Smith MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm); 123404ca555eSLois Curfman McInnes for ( j=0; j<aij->N; j++ ) { 123504ca555eSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 123604ca555eSLois Curfman McInnes } 12370452661fSBarry Smith PetscFree(tmp); PetscFree(tmp2); 123804ca555eSLois Curfman McInnes } 123904ca555eSLois Curfman McInnes else if (type == NORM_INFINITY) { /* max row norm */ 1240515d9167SLois Curfman McInnes double ntemp = 0.0; 124104ca555eSLois Curfman McInnes for ( j=0; j<amat->m; j++ ) { 1242dbb450caSBarry Smith v = amat->a + amat->i[j] + shift; 124304ca555eSLois Curfman McInnes sum = 0.0; 124404ca555eSLois Curfman McInnes for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) { 1245cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 124604ca555eSLois Curfman McInnes } 1247dbb450caSBarry Smith v = bmat->a + bmat->i[j] + shift; 124804ca555eSLois Curfman McInnes for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) { 1249cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 125004ca555eSLois Curfman McInnes } 1251515d9167SLois Curfman McInnes if (sum > ntemp) ntemp = sum; 125204ca555eSLois Curfman McInnes } 12536d84be18SBarry Smith MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm); 125404ca555eSLois Curfman McInnes } 125504ca555eSLois Curfman McInnes else { 1256515d9167SLois Curfman McInnes SETERRQ(1,"MatNorm_MPIAIJ:No support for two norm"); 125704ca555eSLois Curfman McInnes } 125837fa93a5SLois Curfman McInnes } 1259855ac2c5SLois Curfman McInnes return 0; 1260855ac2c5SLois Curfman McInnes } 1261855ac2c5SLois Curfman McInnes 12620de55854SLois Curfman McInnes static int MatTranspose_MPIAIJ(Mat A,Mat *matout) 1263b7c46309SBarry Smith { 1264b7c46309SBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 1265dbb450caSBarry Smith Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data; 1266416022c9SBarry Smith int ierr,shift = Aloc->indexshift; 1267416022c9SBarry Smith Mat B; 1268b7c46309SBarry Smith int M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct; 1269b7c46309SBarry Smith Scalar *array; 1270b7c46309SBarry Smith 12713638b69dSLois Curfman McInnes if (matout == PETSC_NULL && M != N) 12723638b69dSLois Curfman McInnes SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place"); 1273b4fd4287SBarry Smith ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,PETSC_NULL,0, 1274b4fd4287SBarry Smith PETSC_NULL,&B); CHKERRQ(ierr); 1275b7c46309SBarry Smith 1276b7c46309SBarry Smith /* copy over the A part */ 1277ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->A->data; 1278b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1279b7c46309SBarry Smith row = a->rstart; 1280dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;} 1281b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1282416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1283b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 1284b7c46309SBarry Smith } 1285b7c46309SBarry Smith aj = Aloc->j; 12864af08d9eSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;} 1287b7c46309SBarry Smith 1288b7c46309SBarry Smith /* copy over the B part */ 1289ec8511deSBarry Smith Aloc = (Mat_SeqAIJ*) a->B->data; 1290b7c46309SBarry Smith m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a; 1291b7c46309SBarry Smith row = a->rstart; 12920452661fSBarry Smith ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols); 1293dbb450caSBarry Smith for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];} 1294b7c46309SBarry Smith for ( i=0; i<m; i++ ) { 1295416022c9SBarry Smith ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr); 1296b7c46309SBarry Smith row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 1297b7c46309SBarry Smith } 12980452661fSBarry Smith PetscFree(ct); 1299b7c46309SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1300b7c46309SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 13013638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 13020de55854SLois Curfman McInnes *matout = B; 13030de55854SLois Curfman McInnes } else { 13040de55854SLois Curfman McInnes /* This isn't really an in-place transpose .... but free data structures from a */ 13050452661fSBarry Smith PetscFree(a->rowners); 13060de55854SLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 13070de55854SLois Curfman McInnes ierr = MatDestroy(a->B); CHKERRQ(ierr); 13080452661fSBarry Smith if (a->colmap) PetscFree(a->colmap); 13090452661fSBarry Smith if (a->garray) PetscFree(a->garray); 13100de55854SLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 1311a56f8943SBarry Smith if (a->Mvctx) VecScatterDestroy(a->Mvctx); 13120452661fSBarry Smith PetscFree(a); 1313416022c9SBarry Smith PetscMemcpy(A,B,sizeof(struct _Mat)); 13140452661fSBarry Smith PetscHeaderDestroy(B); 13150de55854SLois Curfman McInnes } 1316b7c46309SBarry Smith return 0; 1317b7c46309SBarry Smith } 1318b7c46309SBarry Smith 1319682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat); 1320682d7d0cSBarry Smith static int MatPrintHelp_MPIAIJ(Mat A) 1321682d7d0cSBarry Smith { 1322682d7d0cSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; 1323682d7d0cSBarry Smith 1324682d7d0cSBarry Smith if (!a->rank) return MatPrintHelp_SeqAIJ(a->A); 1325682d7d0cSBarry Smith else return 0; 1326682d7d0cSBarry Smith } 1327682d7d0cSBarry Smith 1328fc76ce87SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MatType,Mat *); 13293d1612f7SBarry Smith static int MatConvertSameType_MPIAIJ(Mat,Mat *,int); 13302f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int); 1331598137ffSSatish Balay int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **); 13328a729477SBarry Smith /* -------------------------------------------------------------------*/ 13332ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ, 133439e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 133544a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 133644a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 1337299609e3SLois Curfman McInnes MatSolve_MPIAIJ,MatSolveAdd_MPIAIJ, 1338299609e3SLois Curfman McInnes MatSolveTrans_MPIAIJ,MatSolveTransAdd_MPIAIJ, 1339299609e3SLois Curfman McInnes MatLUFactor_MPIAIJ,0, 134044a69424SLois Curfman McInnes MatRelax_MPIAIJ, 1341b7c46309SBarry Smith MatTranspose_MPIAIJ, 1342a66be287SLois Curfman McInnes MatGetInfo_MPIAIJ,0, 1343855ac2c5SLois Curfman McInnes MatGetDiagonal_MPIAIJ,0,MatNorm_MPIAIJ, 1344ee50ffe9SBarry Smith MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ, 13451eb62cbbSBarry Smith 0, 1346299609e3SLois Curfman McInnes MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ, 1347299609e3SLois Curfman McInnes MatGetReordering_MPIAIJ, 1348299609e3SLois Curfman McInnes MatLUFactorSymbolic_MPIAIJ,MatLUFactorNumeric_MPIAIJ,0,0, 1349d1710a03SLois Curfman McInnes MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ, 1350299609e3SLois Curfman McInnes MatILUFactorSymbolic_MPIAIJ,0, 13513d1612f7SBarry Smith 0,0,MatConvert_MPIAIJ,0,0,MatConvertSameType_MPIAIJ,0,0, 1352b49de8d1SLois Curfman McInnes 0,0,0, 1353598137ffSSatish Balay MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0, 1354052efed2SBarry Smith MatPrintHelp_MPIAIJ, 1355052efed2SBarry Smith MatScale_MPIAIJ}; 13568a729477SBarry Smith 13571987afe7SBarry Smith /*@C 1358ff756334SLois Curfman McInnes MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format 13593a511b96SLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 13603a511b96SLois Curfman McInnes the user should preallocate the matrix storage by setting the parameters 13613a511b96SLois Curfman McInnes d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately, 13623a511b96SLois Curfman McInnes performance can be increased by more than a factor of 50. 13638a729477SBarry Smith 13648a729477SBarry Smith Input Parameters: 13651eb62cbbSBarry Smith . comm - MPI communicator 13667d3e4905SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 13677d3e4905SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated 13687d3e4905SLois Curfman McInnes if N is given) 13697d3e4905SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 13707d3e4905SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated 13717d3e4905SLois Curfman McInnes if n is given) 1372ab693e5aSLois Curfman McInnes . d_nz - number of nonzeros per row in diagonal portion of local submatrix 1373ff756334SLois Curfman McInnes (same for all local rows) 1374ab693e5aSLois Curfman McInnes . d_nzz - number of nonzeros per row in diagonal portion of local submatrix 1375ab693e5aSLois Curfman McInnes or null (possibly different for each row). You must leave room 1376ab693e5aSLois Curfman McInnes for the diagonal entry even if it is zero. 1377ab693e5aSLois Curfman McInnes . o_nz - number of nonzeros per row in off-diagonal portion of local 1378ab693e5aSLois Curfman McInnes submatrix (same for all local rows). 1379ab693e5aSLois Curfman McInnes . o_nzz - number of nonzeros per row in off-diagonal portion of local 1380ab693e5aSLois Curfman McInnes submatrix or null (possibly different for each row). 13818a729477SBarry Smith 1382ff756334SLois Curfman McInnes Output Parameter: 13838a729477SBarry Smith . newmat - the matrix 13848a729477SBarry Smith 1385ff756334SLois Curfman McInnes Notes: 1386ff756334SLois Curfman McInnes The AIJ format (also called the Yale sparse matrix format or 1387ff756334SLois Curfman McInnes compressed row storage), is fully compatible with standard Fortran 77 13880002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 13890002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 1390ff756334SLois Curfman McInnes 1391ff756334SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 1392ff756334SLois Curfman McInnes (possibly both). 1393ff756334SLois Curfman McInnes 13945511cfe3SLois Curfman McInnes By default, this format uses inodes (identical nodes) when possible. 13955511cfe3SLois Curfman McInnes We search for consecutive rows with the same nonzero structure, thereby 13965511cfe3SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 13975511cfe3SLois Curfman McInnes 13985511cfe3SLois Curfman McInnes Options Database Keys: 13995511cfe3SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 14006e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 14016e62573dSLois Curfman McInnes $ (max limit=5) 14026e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 14036e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 14046e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 14055511cfe3SLois Curfman McInnes 1406e0245417SLois Curfman McInnes Storage Information: 1407e0245417SLois Curfman McInnes For a square global matrix we define each processor's diagonal portion 1408e0245417SLois Curfman McInnes to be its local rows and the corresponding columns (a square submatrix); 1409e0245417SLois Curfman McInnes each processor's off-diagonal portion encompasses the remainder of the 1410e0245417SLois Curfman McInnes local matrix (a rectangular submatrix). 1411e0245417SLois Curfman McInnes 1412e0245417SLois Curfman McInnes The user can specify preallocated storage for the diagonal part of 14135ace5be8SLois Curfman McInnes the local submatrix with either d_nz or d_nnz (not both). Set 14145ace5be8SLois Curfman McInnes d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic 14155ace5be8SLois Curfman McInnes memory allocation. Likewise, specify preallocated storage for the 14165ace5be8SLois Curfman McInnes off-diagonal part of the local submatrix with o_nz or o_nnz (not both). 1417ff756334SLois Curfman McInnes 14185511cfe3SLois Curfman McInnes Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In 14195511cfe3SLois Curfman McInnes the figure below we depict these three local rows and all columns (0-11). 14202191d07cSBarry Smith 1421b810aeb4SBarry Smith $ 0 1 2 3 4 5 6 7 8 9 10 11 1422b810aeb4SBarry Smith $ ------------------- 14235511cfe3SLois Curfman McInnes $ row 3 | o o o d d d o o o o o o 14245511cfe3SLois Curfman McInnes $ row 4 | o o o d d d o o o o o o 14255511cfe3SLois Curfman McInnes $ row 5 | o o o d d d o o o o o o 14265511cfe3SLois Curfman McInnes $ ------------------- 1427b810aeb4SBarry Smith $ 1428b810aeb4SBarry Smith 14295511cfe3SLois Curfman McInnes Thus, any entries in the d locations are stored in the d (diagonal) 14305511cfe3SLois Curfman McInnes submatrix, and any entries in the o locations are stored in the 14315511cfe3SLois Curfman McInnes o (off-diagonal) submatrix. Note that the d and the o submatrices are 14325511cfe3SLois Curfman McInnes stored simply in the MATSEQAIJ format for compressed row storage. 14335511cfe3SLois Curfman McInnes 14345511cfe3SLois Curfman McInnes Now d_nz should indicate the number of nonzeros per row in the d matrix, 14355511cfe3SLois Curfman McInnes and o_nz should indicate the number of nonzeros per row in the o matrix. 14365511cfe3SLois Curfman McInnes In general, for PDE problems in which most nonzeros are near the diagonal, 14373a511b96SLois Curfman McInnes one expects d_nz >> o_nz. For additional details, see the users manual 14383a511b96SLois Curfman McInnes chapter on matrices and the file $(PETSC_DIR)/Performance. 14393a511b96SLois Curfman McInnes 1440dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel 1441ff756334SLois Curfman McInnes 1442fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 14438a729477SBarry Smith @*/ 14441eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 14451eb62cbbSBarry Smith int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *newmat) 14468a729477SBarry Smith { 14478a729477SBarry Smith Mat mat; 1448416022c9SBarry Smith Mat_MPIAIJ *a; 14496abc6512SBarry Smith int ierr, i,sum[2],work[2]; 1450416022c9SBarry Smith 14518a729477SBarry Smith *newmat = 0; 14520452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm); 1453a5a9c739SBarry Smith PLogObjectCreate(mat); 14540452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1455416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 145644a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 145744a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 14588a729477SBarry Smith mat->factor = 0; 1459c456f294SBarry Smith mat->assembled = PETSC_FALSE; 1460d6dfbf8fSBarry Smith 146164119d58SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 146217699dbbSLois Curfman McInnes MPI_Comm_rank(comm,&a->rank); 146317699dbbSLois Curfman McInnes MPI_Comm_size(comm,&a->size); 14641eb62cbbSBarry Smith 1465b4fd4287SBarry Smith if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) 14661987afe7SBarry Smith SETERRQ(1,"MatCreateMPIAIJ:Cannot have PETSc decide rows but set d_nnz or o_nnz"); 14671987afe7SBarry Smith 1468dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE || N == PETSC_DECIDE) { 14691eb62cbbSBarry Smith work[0] = m; work[1] = n; 1470d65a2f8fSBarry Smith MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm ); 1471dbd7a890SLois Curfman McInnes if (M == PETSC_DECIDE) M = sum[0]; 1472dbd7a890SLois Curfman McInnes if (N == PETSC_DECIDE) N = sum[1]; 14731eb62cbbSBarry Smith } 147417699dbbSLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);} 147517699dbbSLois Curfman McInnes if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 1476416022c9SBarry Smith a->m = m; 1477416022c9SBarry Smith a->n = n; 1478416022c9SBarry Smith a->N = N; 1479416022c9SBarry Smith a->M = M; 14801eb62cbbSBarry Smith 14811eb62cbbSBarry Smith /* build local table of row and column ownerships */ 14820452661fSBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1483579c6b6fSBarry Smith PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ)); 148417699dbbSLois Curfman McInnes a->cowners = a->rowners + a->size + 1; 1485416022c9SBarry Smith MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm); 1486416022c9SBarry Smith a->rowners[0] = 0; 148717699dbbSLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 1488416022c9SBarry Smith a->rowners[i] += a->rowners[i-1]; 14898a729477SBarry Smith } 149017699dbbSLois Curfman McInnes a->rstart = a->rowners[a->rank]; 149117699dbbSLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 1492416022c9SBarry Smith MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm); 1493416022c9SBarry Smith a->cowners[0] = 0; 149417699dbbSLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 1495416022c9SBarry Smith a->cowners[i] += a->cowners[i-1]; 14968a729477SBarry Smith } 149717699dbbSLois Curfman McInnes a->cstart = a->cowners[a->rank]; 149817699dbbSLois Curfman McInnes a->cend = a->cowners[a->rank+1]; 14998a729477SBarry Smith 15005ace5be8SLois Curfman McInnes if (d_nz == PETSC_DEFAULT) d_nz = 5; 1501416022c9SBarry Smith ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&a->A); CHKERRQ(ierr); 1502416022c9SBarry Smith PLogObjectParent(mat,a->A); 15037b8455f0SLois Curfman McInnes if (o_nz == PETSC_DEFAULT) o_nz = 0; 1504416022c9SBarry Smith ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&a->B); CHKERRQ(ierr); 1505416022c9SBarry Smith PLogObjectParent(mat,a->B); 15068a729477SBarry Smith 15071eb62cbbSBarry Smith /* build cache for off array entries formed */ 1508416022c9SBarry Smith ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr); 1509416022c9SBarry Smith a->colmap = 0; 1510416022c9SBarry Smith a->garray = 0; 15114b0e389bSBarry Smith a->roworiented = 1; 15128a729477SBarry Smith 15131eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 1514416022c9SBarry Smith a->lvec = 0; 1515416022c9SBarry Smith a->Mvctx = 0; 15168a729477SBarry Smith 15177a0afa10SBarry Smith /* stuff for MatGetRow() */ 15187a0afa10SBarry Smith a->rowindices = 0; 15197a0afa10SBarry Smith a->rowvalues = 0; 15207a0afa10SBarry Smith a->getrowactive = PETSC_FALSE; 15217a0afa10SBarry Smith 1522d6dfbf8fSBarry Smith *newmat = mat; 1523d6dfbf8fSBarry Smith return 0; 1524d6dfbf8fSBarry Smith } 1525c74985f6SBarry Smith 15263d1612f7SBarry Smith static int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues) 1527d6dfbf8fSBarry Smith { 1528d6dfbf8fSBarry Smith Mat mat; 1529416022c9SBarry Smith Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data; 153025cdf11fSBarry Smith int ierr, len,flg; 1531d6dfbf8fSBarry Smith 1532416022c9SBarry Smith *newmat = 0; 15330452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm); 1534a5a9c739SBarry Smith PLogObjectCreate(mat); 15350452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a); 1536416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 153744a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 153844a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1539d6dfbf8fSBarry Smith mat->factor = matin->factor; 1540c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1541d6dfbf8fSBarry Smith 1542416022c9SBarry Smith a->m = oldmat->m; 1543416022c9SBarry Smith a->n = oldmat->n; 1544416022c9SBarry Smith a->M = oldmat->M; 1545416022c9SBarry Smith a->N = oldmat->N; 1546d6dfbf8fSBarry Smith 1547416022c9SBarry Smith a->rstart = oldmat->rstart; 1548416022c9SBarry Smith a->rend = oldmat->rend; 1549416022c9SBarry Smith a->cstart = oldmat->cstart; 1550416022c9SBarry Smith a->cend = oldmat->cend; 155117699dbbSLois Curfman McInnes a->size = oldmat->size; 155217699dbbSLois Curfman McInnes a->rank = oldmat->rank; 155364119d58SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 1554*bcd2baecSBarry Smith a->rowvalues = 0; 1555*bcd2baecSBarry Smith a->getrowactive = PETSC_FALSE; 1556d6dfbf8fSBarry Smith 15570452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 155817699dbbSLois Curfman McInnes PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ)); 155917699dbbSLois Curfman McInnes PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 1560416022c9SBarry Smith ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 15612ee70a88SLois Curfman McInnes if (oldmat->colmap) { 15620452661fSBarry Smith a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap); 1563416022c9SBarry Smith PLogObjectMemory(mat,(a->N)*sizeof(int)); 1564416022c9SBarry Smith PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int)); 1565416022c9SBarry Smith } else a->colmap = 0; 1566ec8511deSBarry Smith if (oldmat->garray && (len = ((Mat_SeqAIJ *) (oldmat->B->data))->n)) { 15670452661fSBarry Smith a->garray = (int *) PetscMalloc(len*sizeof(int)); CHKPTRQ(a->garray); 1568464493b3SBarry Smith PLogObjectMemory(mat,len*sizeof(int)); 1569416022c9SBarry Smith PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int)); 1570416022c9SBarry Smith } else a->garray = 0; 1571d6dfbf8fSBarry Smith 1572416022c9SBarry Smith ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1573416022c9SBarry Smith PLogObjectParent(mat,a->lvec); 1574a56f8943SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1575416022c9SBarry Smith PLogObjectParent(mat,a->Mvctx); 1576416022c9SBarry Smith ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 1577416022c9SBarry Smith PLogObjectParent(mat,a->A); 1578416022c9SBarry Smith ierr = MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr); 1579416022c9SBarry Smith PLogObjectParent(mat,a->B); 15805dd7a6c7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 158125cdf11fSBarry Smith if (flg) { 1582682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1583682d7d0cSBarry Smith } 15848a729477SBarry Smith *newmat = mat; 15858a729477SBarry Smith return 0; 15868a729477SBarry Smith } 1587416022c9SBarry Smith 1588416022c9SBarry Smith #include "sysio.h" 1589416022c9SBarry Smith 1590c456f294SBarry Smith int MatLoad_MPIAIJ(Viewer bview,MatType type,Mat *newmat) 1591416022c9SBarry Smith { 1592d65a2f8fSBarry Smith Mat A; 1593416022c9SBarry Smith int i, nz, ierr, j,rstart, rend, fd; 1594d65a2f8fSBarry Smith Scalar *vals,*svals; 1595416022c9SBarry Smith PetscObject vobj = (PetscObject) bview; 1596416022c9SBarry Smith MPI_Comm comm = vobj->comm; 1597416022c9SBarry Smith MPI_Status status; 159817699dbbSLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1599d65a2f8fSBarry Smith int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 1600d65a2f8fSBarry Smith int tag = ((PetscObject)bview)->tag; 1601416022c9SBarry Smith 160217699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 160317699dbbSLois Curfman McInnes if (!rank) { 1604416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1605416022c9SBarry Smith ierr = SYRead(fd,(char *)header,4,SYINT); CHKERRQ(ierr); 1606c456f294SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIAIJ:not matrix object"); 1607416022c9SBarry Smith } 1608416022c9SBarry Smith 1609416022c9SBarry Smith MPI_Bcast(header+1,3,MPI_INT,0,comm); 1610416022c9SBarry Smith M = header[1]; N = header[2]; 1611416022c9SBarry Smith /* determine ownership of all rows */ 161217699dbbSLois Curfman McInnes m = M/size + ((M % size) > rank); 16130452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1614416022c9SBarry Smith MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 1615416022c9SBarry Smith rowners[0] = 0; 161617699dbbSLois Curfman McInnes for ( i=2; i<=size; i++ ) { 1617416022c9SBarry Smith rowners[i] += rowners[i-1]; 1618416022c9SBarry Smith } 161917699dbbSLois Curfman McInnes rstart = rowners[rank]; 162017699dbbSLois Curfman McInnes rend = rowners[rank+1]; 1621416022c9SBarry Smith 1622416022c9SBarry Smith /* distribute row lengths to all processors */ 16230452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1624416022c9SBarry Smith offlens = ourlens + (rend-rstart); 162517699dbbSLois Curfman McInnes if (!rank) { 16260452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1627416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 16280452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 162917699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1630d65a2f8fSBarry Smith MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 16310452661fSBarry Smith PetscFree(sndcounts); 1632416022c9SBarry Smith } 1633416022c9SBarry Smith else { 1634416022c9SBarry Smith MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 1635416022c9SBarry Smith } 1636416022c9SBarry Smith 163717699dbbSLois Curfman McInnes if (!rank) { 1638416022c9SBarry Smith /* calculate the number of nonzeros on each processor */ 16390452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1640cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 164117699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 1642416022c9SBarry Smith for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1643416022c9SBarry Smith procsnz[i] += rowlengths[j]; 1644416022c9SBarry Smith } 1645416022c9SBarry Smith } 16460452661fSBarry Smith PetscFree(rowlengths); 1647416022c9SBarry Smith 1648416022c9SBarry Smith /* determine max buffer needed and allocate it */ 1649416022c9SBarry Smith maxnz = 0; 165017699dbbSLois Curfman McInnes for ( i=0; i<size; i++ ) { 16510452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 1652416022c9SBarry Smith } 16530452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1654416022c9SBarry Smith 1655416022c9SBarry Smith /* read in my part of the matrix column indices */ 1656416022c9SBarry Smith nz = procsnz[0]; 16570452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1658d65a2f8fSBarry Smith ierr = SYRead(fd,mycols,nz,SYINT); CHKERRQ(ierr); 1659d65a2f8fSBarry Smith 1660d65a2f8fSBarry Smith /* read in every one elses and ship off */ 166117699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1662d65a2f8fSBarry Smith nz = procsnz[i]; 1663416022c9SBarry Smith ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr); 1664d65a2f8fSBarry Smith MPI_Send(cols,nz,MPI_INT,i,tag,comm); 1665d65a2f8fSBarry Smith } 16660452661fSBarry Smith PetscFree(cols); 1667416022c9SBarry Smith } 1668416022c9SBarry Smith else { 1669416022c9SBarry Smith /* determine buffer space needed for message */ 1670416022c9SBarry Smith nz = 0; 1671416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1672416022c9SBarry Smith nz += ourlens[i]; 1673416022c9SBarry Smith } 16740452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1675416022c9SBarry Smith 1676416022c9SBarry Smith /* receive message of column indices*/ 1677d65a2f8fSBarry Smith MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 1678416022c9SBarry Smith MPI_Get_count(&status,MPI_INT,&maxnz); 1679c456f294SBarry Smith if (maxnz != nz) SETERRQ(1,"MatLoad_MPIAIJ:something is wrong with file"); 1680416022c9SBarry Smith } 1681416022c9SBarry Smith 1682416022c9SBarry Smith /* loop over local rows, determining number of off diagonal entries */ 1683cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 1684416022c9SBarry Smith jj = 0; 1685416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1686416022c9SBarry Smith for ( j=0; j<ourlens[i]; j++ ) { 1687d65a2f8fSBarry Smith if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1688416022c9SBarry Smith jj++; 1689416022c9SBarry Smith } 1690416022c9SBarry Smith } 1691d65a2f8fSBarry Smith 1692d65a2f8fSBarry Smith /* create our matrix */ 1693416022c9SBarry Smith for ( i=0; i<m; i++ ) { 1694416022c9SBarry Smith ourlens[i] -= offlens[i]; 1695416022c9SBarry Smith } 1696d65a2f8fSBarry Smith ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr); 1697d65a2f8fSBarry Smith A = *newmat; 1698d65a2f8fSBarry Smith MatSetOption(A,COLUMNS_SORTED); 1699d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1700d65a2f8fSBarry Smith ourlens[i] += offlens[i]; 1701d65a2f8fSBarry Smith } 1702416022c9SBarry Smith 170317699dbbSLois Curfman McInnes if (!rank) { 17040452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1705416022c9SBarry Smith 1706416022c9SBarry Smith /* read in my part of the matrix numerical values */ 1707416022c9SBarry Smith nz = procsnz[0]; 1708416022c9SBarry Smith ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 1709d65a2f8fSBarry Smith 1710d65a2f8fSBarry Smith /* insert into matrix */ 1711d65a2f8fSBarry Smith jj = rstart; 1712d65a2f8fSBarry Smith smycols = mycols; 1713d65a2f8fSBarry Smith svals = vals; 1714d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1715d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1716d65a2f8fSBarry Smith smycols += ourlens[i]; 1717d65a2f8fSBarry Smith svals += ourlens[i]; 1718d65a2f8fSBarry Smith jj++; 1719416022c9SBarry Smith } 1720416022c9SBarry Smith 1721d65a2f8fSBarry Smith /* read in other processors and ship out */ 172217699dbbSLois Curfman McInnes for ( i=1; i<size; i++ ) { 1723416022c9SBarry Smith nz = procsnz[i]; 1724416022c9SBarry Smith ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 1725d65a2f8fSBarry Smith MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 1726416022c9SBarry Smith } 17270452661fSBarry Smith PetscFree(procsnz); 1728416022c9SBarry Smith } 1729d65a2f8fSBarry Smith else { 1730d65a2f8fSBarry Smith /* receive numeric values */ 17310452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1732416022c9SBarry Smith 1733d65a2f8fSBarry Smith /* receive message of values*/ 1734d65a2f8fSBarry Smith MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 1735d65a2f8fSBarry Smith MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 1736c456f294SBarry Smith if (maxnz != nz) SETERRQ(1,"MatLoad_MPIAIJ:something is wrong with file"); 1737d65a2f8fSBarry Smith 1738d65a2f8fSBarry Smith /* insert into matrix */ 1739d65a2f8fSBarry Smith jj = rstart; 1740d65a2f8fSBarry Smith smycols = mycols; 1741d65a2f8fSBarry Smith svals = vals; 1742d65a2f8fSBarry Smith for ( i=0; i<m; i++ ) { 1743d65a2f8fSBarry Smith ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1744d65a2f8fSBarry Smith smycols += ourlens[i]; 1745d65a2f8fSBarry Smith svals += ourlens[i]; 1746d65a2f8fSBarry Smith jj++; 1747d65a2f8fSBarry Smith } 1748d65a2f8fSBarry Smith } 17490452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1750d65a2f8fSBarry Smith 1751d65a2f8fSBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 1752d65a2f8fSBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 1753416022c9SBarry Smith return 0; 1754416022c9SBarry Smith } 1755