xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision 4d39ef84ea3d244593663a951516917320ccfd4b)
1cb512458SBarry Smith #ifndef lint
2*4d39ef84SLois Curfman McInnes static char vcid[] = "$Id: mpiaij.c,v 1.88 1995/10/12 13:41:03 curfman Exp curfman $";
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 
2078b31e54SBarry Smith   aij->colmap = (int *) PETSCMALLOC(aij->N*sizeof(int));CHKPTRQ(aij->colmap);
21464493b3SBarry Smith   PLogObjectMemory(mat,aij->N*sizeof(int));
22416022c9SBarry Smith   PetscZero(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;
33299609e3SLois Curfman McInnes   if (aij->numtids == 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 
392ee70a88SLois Curfman McInnes static int MatSetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,
401eb62cbbSBarry Smith                                int *idxn,Scalar *v,InsertMode addv)
418a729477SBarry Smith {
4244a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
43dbb450caSBarry Smith   Mat_SeqAIJ *C = (Mat_SeqAIJ*) aij->A->data;
441eb62cbbSBarry Smith   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
451eb62cbbSBarry Smith   int        cstart = aij->cstart, cend = aij->cend,row,col;
46dbb450caSBarry Smith   int        shift = C->indexshift;
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++ ) {
53bbb6d6a8SBarry Smith     if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIAIJ:Negative row");
54bbb6d6a8SBarry Smith     if (idxm[i] >= aij->M) SETERRQ(1,"MatSetValues_MPIAIJ:Row too large");
551eb62cbbSBarry Smith     if (idxm[i] >= rstart && idxm[i] < rend) {
561eb62cbbSBarry Smith       row = idxm[i] - rstart;
571eb62cbbSBarry Smith       for ( j=0; j<n; j++ ) {
58bbb6d6a8SBarry Smith         if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIAIJ:Negative column");
59bbb6d6a8SBarry Smith         if (idxn[j] >= aij->N) SETERRQ(1,"MatSetValues_MPIAIJ:Col too large");
601eb62cbbSBarry Smith         if (idxn[j] >= cstart && idxn[j] < cend){
611eb62cbbSBarry Smith           col = idxn[j] - cstart;
6278b31e54SBarry Smith           ierr = MatSetValues(aij->A,1,&row,1,&col,v+i*n+j,addv);CHKERRQ(ierr);
631eb62cbbSBarry Smith         }
641eb62cbbSBarry Smith         else {
65d6dfbf8fSBarry Smith           if (aij->assembled) {
66b7c46309SBarry Smith             if (!aij->colmap) {ierr = CreateColmap_Private(mat);CHKERRQ(ierr);}
67dbb450caSBarry Smith             col = aij->colmap[idxn[j]] + shift;
68ec8511deSBarry Smith             if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
692493cbb0SBarry Smith               ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
702493cbb0SBarry Smith               col =  idxn[j];
71d6dfbf8fSBarry Smith             }
729e25ed09SBarry Smith           }
739e25ed09SBarry Smith           else col = idxn[j];
7478b31e54SBarry Smith           ierr = MatSetValues(aij->B,1,&row,1,&col,v+i*n+j,addv);CHKERRQ(ierr);
751eb62cbbSBarry Smith         }
761eb62cbbSBarry Smith       }
771eb62cbbSBarry Smith     }
781eb62cbbSBarry Smith     else {
79416022c9SBarry Smith       ierr = StashValues_Private(&aij->stash,idxm[i],n,idxn,v+i*n,addv);CHKERRQ(ierr);
801eb62cbbSBarry Smith     }
818a729477SBarry Smith   }
828a729477SBarry Smith   return 0;
838a729477SBarry Smith }
848a729477SBarry Smith 
85fc76ce87SLois Curfman McInnes static int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode)
868a729477SBarry Smith {
8744a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
88d6dfbf8fSBarry Smith   MPI_Comm    comm = mat->comm;
896abc6512SBarry Smith   int         numtids = aij->numtids, *owners = aij->rowners;
90416022c9SBarry Smith   int         mytid = aij->mytid,tag = mat->tag, *owner,*starts,count,ierr;
911eb62cbbSBarry Smith   MPI_Request *send_waits,*recv_waits;
926abc6512SBarry Smith   int         *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
931eb62cbbSBarry Smith   InsertMode  addv;
941eb62cbbSBarry Smith   Scalar      *rvalues,*svalues;
951eb62cbbSBarry Smith 
961eb62cbbSBarry Smith   /* make sure all processors are either in INSERTMODE or ADDMODE */
97d65a2f8fSBarry Smith   MPI_Allreduce(&aij->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
98dbb450caSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
99bbb6d6a8SBarry Smith     SETERRQ(1,"MatAssemblyBegin_MPIAIJ:Some processors inserted others added");
1001eb62cbbSBarry Smith   }
1011eb62cbbSBarry Smith   aij->insertmode = addv; /* in case this processor had no cache */
1021eb62cbbSBarry Smith 
1031eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
10478b31e54SBarry Smith   nprocs = (int *) PETSCMALLOC( 2*numtids*sizeof(int) ); CHKPTRQ(nprocs);
105416022c9SBarry Smith   PetscZero(nprocs,2*numtids*sizeof(int)); procs = nprocs + numtids;
10678b31e54SBarry Smith   owner = (int *) PETSCMALLOC( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
1071eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
1081eb62cbbSBarry Smith     idx = aij->stash.idx[i];
1091eb62cbbSBarry Smith     for ( j=0; j<numtids; j++ ) {
1101eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
1111eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
1128a729477SBarry Smith       }
1138a729477SBarry Smith     }
1148a729477SBarry Smith   }
1151eb62cbbSBarry Smith   nsends = 0;  for ( i=0; i<numtids; i++ ) { nsends += procs[i];}
1161eb62cbbSBarry Smith 
1171eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
11878b31e54SBarry Smith   work = (int *) PETSCMALLOC( numtids*sizeof(int) ); CHKPTRQ(work);
119d65a2f8fSBarry Smith   MPI_Allreduce(procs, work,numtids,MPI_INT,MPI_SUM,comm);
1201eb62cbbSBarry Smith   nreceives = work[mytid];
121d65a2f8fSBarry Smith   MPI_Allreduce( nprocs, work,numtids,MPI_INT,MPI_MAX,comm);
1221eb62cbbSBarry Smith   nmax = work[mytid];
12378b31e54SBarry Smith   PETSCFREE(work);
1241eb62cbbSBarry Smith 
1251eb62cbbSBarry Smith   /* post receives:
1261eb62cbbSBarry Smith        1) each message will consist of ordered pairs
1271eb62cbbSBarry Smith      (global index,value) we store the global index as a double
128d6dfbf8fSBarry Smith      to simplify the message passing.
1291eb62cbbSBarry Smith        2) since we don't know how long each individual message is we
1301eb62cbbSBarry Smith      allocate the largest needed buffer for each receive. Potentially
1311eb62cbbSBarry Smith      this is a lot of wasted space.
1321eb62cbbSBarry Smith 
1331eb62cbbSBarry Smith 
1341eb62cbbSBarry Smith        This could be done better.
1351eb62cbbSBarry Smith   */
13678b31e54SBarry Smith   rvalues = (Scalar *) PETSCMALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));
13778b31e54SBarry Smith   CHKPTRQ(rvalues);
13878b31e54SBarry Smith   recv_waits = (MPI_Request *) PETSCMALLOC((nreceives+1)*sizeof(MPI_Request));
13978b31e54SBarry Smith   CHKPTRQ(recv_waits);
1401eb62cbbSBarry Smith   for ( i=0; i<nreceives; i++ ) {
141416022c9SBarry Smith     MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,
1421eb62cbbSBarry Smith               comm,recv_waits+i);
1431eb62cbbSBarry Smith   }
1441eb62cbbSBarry Smith 
1451eb62cbbSBarry Smith   /* do sends:
1461eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
1471eb62cbbSBarry Smith          the ith processor
1481eb62cbbSBarry Smith   */
149416022c9SBarry Smith   svalues = (Scalar *) PETSCMALLOC(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues);
15078b31e54SBarry Smith   send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request));
15178b31e54SBarry Smith   CHKPTRQ(send_waits);
15278b31e54SBarry Smith   starts = (int *) PETSCMALLOC( numtids*sizeof(int) ); CHKPTRQ(starts);
1531eb62cbbSBarry Smith   starts[0] = 0;
1541eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1551eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
1561eb62cbbSBarry Smith     svalues[3*starts[owner[i]]]       = (Scalar)  aij->stash.idx[i];
1571eb62cbbSBarry Smith     svalues[3*starts[owner[i]]+1]     = (Scalar)  aij->stash.idy[i];
1581eb62cbbSBarry Smith     svalues[3*(starts[owner[i]]++)+2] =  aij->stash.array[i];
1591eb62cbbSBarry Smith   }
16078b31e54SBarry Smith   PETSCFREE(owner);
1611eb62cbbSBarry Smith   starts[0] = 0;
1621eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1631eb62cbbSBarry Smith   count = 0;
1641eb62cbbSBarry Smith   for ( i=0; i<numtids; i++ ) {
1651eb62cbbSBarry Smith     if (procs[i]) {
166416022c9SBarry Smith       MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,
1671eb62cbbSBarry Smith                 comm,send_waits+count++);
1681eb62cbbSBarry Smith     }
1691eb62cbbSBarry Smith   }
17078b31e54SBarry Smith   PETSCFREE(starts); PETSCFREE(nprocs);
1711eb62cbbSBarry Smith 
1721eb62cbbSBarry Smith   /* Free cache space */
17378b31e54SBarry Smith   ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr);
1741eb62cbbSBarry Smith 
1751eb62cbbSBarry Smith   aij->svalues    = svalues;    aij->rvalues    = rvalues;
1761eb62cbbSBarry Smith   aij->nsends     = nsends;     aij->nrecvs     = nreceives;
1771eb62cbbSBarry Smith   aij->send_waits = send_waits; aij->recv_waits = recv_waits;
1781eb62cbbSBarry Smith   aij->rmax       = nmax;
1791eb62cbbSBarry Smith 
1801eb62cbbSBarry Smith   return 0;
1811eb62cbbSBarry Smith }
18244a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat);
1831eb62cbbSBarry Smith 
184fc76ce87SLois Curfman McInnes static int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode)
1851eb62cbbSBarry Smith {
18644a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
187dbb450caSBarry Smith   Mat_SeqAIJ *C = (Mat_SeqAIJ *) aij->A->data;
1881eb62cbbSBarry Smith   MPI_Status  *send_status,recv_status;
189416022c9SBarry Smith   int         imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr;
190416022c9SBarry Smith   int         row,col,other_disassembled,shift = C->indexshift;
1911eb62cbbSBarry Smith   Scalar      *values,val;
1921eb62cbbSBarry Smith   InsertMode  addv = aij->insertmode;
1931eb62cbbSBarry Smith 
1941eb62cbbSBarry Smith   /*  wait on receives */
1951eb62cbbSBarry Smith   while (count) {
196d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);
1971eb62cbbSBarry Smith     /* unpack receives into our local space */
198d6dfbf8fSBarry Smith     values = aij->rvalues + 3*imdex*aij->rmax;
199c60448a5SBarry Smith     MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
2001eb62cbbSBarry Smith     n = n/3;
2011eb62cbbSBarry Smith     for ( i=0; i<n; i++ ) {
2021eb62cbbSBarry Smith       row = (int) PETSCREAL(values[3*i]) - aij->rstart;
2031eb62cbbSBarry Smith       col = (int) PETSCREAL(values[3*i+1]);
2041eb62cbbSBarry Smith       val = values[3*i+2];
2051eb62cbbSBarry Smith       if (col >= aij->cstart && col < aij->cend) {
2061eb62cbbSBarry Smith         col -= aij->cstart;
2071eb62cbbSBarry Smith         MatSetValues(aij->A,1,&row,1,&col,&val,addv);
2081eb62cbbSBarry Smith       }
2091eb62cbbSBarry Smith       else {
210d6dfbf8fSBarry Smith         if (aij->assembled) {
211b7c46309SBarry Smith           if (!aij->colmap) {ierr = CreateColmap_Private(mat);CHKERRQ(ierr);}
212dbb450caSBarry Smith           col = aij->colmap[col] + shift;
213ec8511deSBarry Smith           if (col < 0  && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
2142493cbb0SBarry Smith             ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
2152493cbb0SBarry Smith             col = (int) PETSCREAL(values[3*i+1]);
216d6dfbf8fSBarry Smith           }
2179e25ed09SBarry Smith         }
2181eb62cbbSBarry Smith         MatSetValues(aij->B,1,&row,1,&col,&val,addv);
2191eb62cbbSBarry Smith       }
2201eb62cbbSBarry Smith     }
2211eb62cbbSBarry Smith     count--;
2221eb62cbbSBarry Smith   }
22378b31e54SBarry Smith   PETSCFREE(aij->recv_waits); PETSCFREE(aij->rvalues);
2241eb62cbbSBarry Smith 
2251eb62cbbSBarry Smith   /* wait on sends */
2261eb62cbbSBarry Smith   if (aij->nsends) {
22778b31e54SBarry Smith     send_status = (MPI_Status *) PETSCMALLOC(aij->nsends*sizeof(MPI_Status));
22878b31e54SBarry Smith     CHKPTRQ(send_status);
2291eb62cbbSBarry Smith     MPI_Waitall(aij->nsends,aij->send_waits,send_status);
23078b31e54SBarry Smith     PETSCFREE(send_status);
2311eb62cbbSBarry Smith   }
23278b31e54SBarry Smith   PETSCFREE(aij->send_waits); PETSCFREE(aij->svalues);
2331eb62cbbSBarry Smith 
23464119d58SLois Curfman McInnes   aij->insertmode = NOT_SET_VALUES;
23578b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr);
23678b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr);
2371eb62cbbSBarry Smith 
2382493cbb0SBarry Smith   /* determine if any processor has disassembled, if so we must
2392493cbb0SBarry Smith      also disassemble ourselfs, in order that we may reassemble. */
240416022c9SBarry Smith   MPI_Allreduce(&aij->assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);
2412493cbb0SBarry Smith   if (aij->assembled && !other_disassembled) {
2422493cbb0SBarry Smith     ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
2432493cbb0SBarry Smith   }
2442493cbb0SBarry Smith 
2455e42470aSBarry Smith   if (!aij->assembled && mode == FINAL_ASSEMBLY) {
24678b31e54SBarry Smith     ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr);
247d6dfbf8fSBarry Smith     aij->assembled = 1;
2485e42470aSBarry Smith   }
24978b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr);
25078b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr);
2515e42470aSBarry Smith 
2528a729477SBarry Smith   return 0;
2538a729477SBarry Smith }
2548a729477SBarry Smith 
2552ee70a88SLois Curfman McInnes static int MatZeroEntries_MPIAIJ(Mat A)
2561eb62cbbSBarry Smith {
25744a69424SLois Curfman McInnes   Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data;
258dbd7a890SLois Curfman McInnes   int        ierr;
25978b31e54SBarry Smith   ierr = MatZeroEntries(l->A); CHKERRQ(ierr);
26078b31e54SBarry Smith   ierr = MatZeroEntries(l->B); CHKERRQ(ierr);
2611eb62cbbSBarry Smith   return 0;
2621eb62cbbSBarry Smith }
2631eb62cbbSBarry Smith 
2641eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the
2651eb62cbbSBarry Smith    matrix is square and the column and row owerships are identical.
2661eb62cbbSBarry Smith    This is a BUG. The only way to fix it seems to be to access
2671eb62cbbSBarry Smith    aij->A and aij->B directly and not through the MatZeroRows()
2681eb62cbbSBarry Smith    routine.
2691eb62cbbSBarry Smith */
27044a69424SLois Curfman McInnes static int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag)
2711eb62cbbSBarry Smith {
27244a69424SLois Curfman McInnes   Mat_MPIAIJ     *l = (Mat_MPIAIJ *) A->data;
2731eb62cbbSBarry Smith   int            i,ierr,N, *rows,*owners = l->rowners,numtids = l->numtids;
2746abc6512SBarry Smith   int            *procs,*nprocs,j,found,idx,nsends,*work;
2751eb62cbbSBarry Smith   int            nmax,*svalues,*starts,*owner,nrecvs,mytid = l->mytid;
2765392566eSBarry Smith   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
277d6dfbf8fSBarry Smith   int            *lens,imdex,*lrows,*values;
278d6dfbf8fSBarry Smith   MPI_Comm       comm = A->comm;
2791eb62cbbSBarry Smith   MPI_Request    *send_waits,*recv_waits;
2801eb62cbbSBarry Smith   MPI_Status     recv_status,*send_status;
2811eb62cbbSBarry Smith   IS             istmp;
2821eb62cbbSBarry Smith 
28348d91487SBarry Smith   if (!l->assembled) SETERRQ(1,"MatZeroRows_MPIAIJ:Must assemble matrix");
28478b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
28578b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
2861eb62cbbSBarry Smith 
2871eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
28878b31e54SBarry Smith   nprocs = (int *) PETSCMALLOC( 2*numtids*sizeof(int) ); CHKPTRQ(nprocs);
289416022c9SBarry Smith   PetscZero(nprocs,2*numtids*sizeof(int)); procs = nprocs + numtids;
29078b31e54SBarry Smith   owner = (int *) PETSCMALLOC((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
2911eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
2921eb62cbbSBarry Smith     idx = rows[i];
2931eb62cbbSBarry Smith     found = 0;
2941eb62cbbSBarry Smith     for ( j=0; j<numtids; j++ ) {
2951eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
2961eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
2971eb62cbbSBarry Smith       }
2981eb62cbbSBarry Smith     }
299bbb6d6a8SBarry Smith     if (!found) SETERRQ(1,"MatZeroRows_MPIAIJ:Index out of range");
3001eb62cbbSBarry Smith   }
3011eb62cbbSBarry Smith   nsends = 0;  for ( i=0; i<numtids; i++ ) { nsends += procs[i];}
3021eb62cbbSBarry Smith 
3031eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
30478b31e54SBarry Smith   work = (int *) PETSCMALLOC( numtids*sizeof(int) ); CHKPTRQ(work);
305d65a2f8fSBarry Smith   MPI_Allreduce( procs, work,numtids,MPI_INT,MPI_SUM,comm);
3061eb62cbbSBarry Smith   nrecvs = work[mytid];
307d65a2f8fSBarry Smith   MPI_Allreduce( nprocs, work,numtids,MPI_INT,MPI_MAX,comm);
3081eb62cbbSBarry Smith   nmax = work[mytid];
30978b31e54SBarry Smith   PETSCFREE(work);
3101eb62cbbSBarry Smith 
3111eb62cbbSBarry Smith   /* post receives:   */
31278b31e54SBarry Smith   rvalues = (int *) PETSCMALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
31378b31e54SBarry Smith   CHKPTRQ(rvalues);
31478b31e54SBarry Smith   recv_waits = (MPI_Request *) PETSCMALLOC((nrecvs+1)*sizeof(MPI_Request));
31578b31e54SBarry Smith   CHKPTRQ(recv_waits);
3161eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
317416022c9SBarry Smith     MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
3181eb62cbbSBarry Smith   }
3191eb62cbbSBarry Smith 
3201eb62cbbSBarry Smith   /* do sends:
3211eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
3221eb62cbbSBarry Smith          the ith processor
3231eb62cbbSBarry Smith   */
32478b31e54SBarry Smith   svalues = (int *) PETSCMALLOC( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
32578b31e54SBarry Smith   send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request));
32678b31e54SBarry Smith   CHKPTRQ(send_waits);
32778b31e54SBarry Smith   starts = (int *) PETSCMALLOC( (numtids+1)*sizeof(int) ); CHKPTRQ(starts);
3281eb62cbbSBarry Smith   starts[0] = 0;
3291eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3301eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
3311eb62cbbSBarry Smith     svalues[starts[owner[i]]++] = rows[i];
3321eb62cbbSBarry Smith   }
3331eb62cbbSBarry Smith   ISRestoreIndices(is,&rows);
3341eb62cbbSBarry Smith 
3351eb62cbbSBarry Smith   starts[0] = 0;
3361eb62cbbSBarry Smith   for ( i=1; i<numtids+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3371eb62cbbSBarry Smith   count = 0;
3381eb62cbbSBarry Smith   for ( i=0; i<numtids; i++ ) {
3391eb62cbbSBarry Smith     if (procs[i]) {
340416022c9SBarry Smith       MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);
3411eb62cbbSBarry Smith     }
3421eb62cbbSBarry Smith   }
34378b31e54SBarry Smith   PETSCFREE(starts);
3441eb62cbbSBarry Smith 
3451eb62cbbSBarry Smith   base = owners[mytid];
3461eb62cbbSBarry Smith 
3471eb62cbbSBarry Smith   /*  wait on receives */
34878b31e54SBarry Smith   lens   = (int *) PETSCMALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
3491eb62cbbSBarry Smith   source = lens + nrecvs;
3501eb62cbbSBarry Smith   count  = nrecvs; slen = 0;
3511eb62cbbSBarry Smith   while (count) {
352d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
3531eb62cbbSBarry Smith     /* unpack receives into our local space */
3541eb62cbbSBarry Smith     MPI_Get_count(&recv_status,MPI_INT,&n);
355d6dfbf8fSBarry Smith     source[imdex]  = recv_status.MPI_SOURCE;
356d6dfbf8fSBarry Smith     lens[imdex]  = n;
3571eb62cbbSBarry Smith     slen += n;
3581eb62cbbSBarry Smith     count--;
3591eb62cbbSBarry Smith   }
36078b31e54SBarry Smith   PETSCFREE(recv_waits);
3611eb62cbbSBarry Smith 
3621eb62cbbSBarry Smith   /* move the data into the send scatter */
36378b31e54SBarry Smith   lrows = (int *) PETSCMALLOC( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
3641eb62cbbSBarry Smith   count = 0;
3651eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
3661eb62cbbSBarry Smith     values = rvalues + i*nmax;
3671eb62cbbSBarry Smith     for ( j=0; j<lens[i]; j++ ) {
3681eb62cbbSBarry Smith       lrows[count++] = values[j] - base;
3691eb62cbbSBarry Smith     }
3701eb62cbbSBarry Smith   }
37178b31e54SBarry Smith   PETSCFREE(rvalues); PETSCFREE(lens);
37278b31e54SBarry Smith   PETSCFREE(owner); PETSCFREE(nprocs);
3731eb62cbbSBarry Smith 
3741eb62cbbSBarry Smith   /* actually zap the local rows */
375416022c9SBarry Smith   ierr = ISCreateSeq(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
376464493b3SBarry Smith   PLogObjectParent(A,istmp);
377416022c9SBarry Smith   PETSCFREE(lrows);
37878b31e54SBarry Smith   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
37978b31e54SBarry Smith   ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr);
38078b31e54SBarry Smith   ierr = ISDestroy(istmp); CHKERRQ(ierr);
3811eb62cbbSBarry Smith 
3821eb62cbbSBarry Smith   /* wait on sends */
3831eb62cbbSBarry Smith   if (nsends) {
38478b31e54SBarry Smith     send_status = (MPI_Status *) PETSCMALLOC(nsends*sizeof(MPI_Status));
38578b31e54SBarry Smith     CHKPTRQ(send_status);
3861eb62cbbSBarry Smith     MPI_Waitall(nsends,send_waits,send_status);
38778b31e54SBarry Smith     PETSCFREE(send_status);
3881eb62cbbSBarry Smith   }
38978b31e54SBarry Smith   PETSCFREE(send_waits); PETSCFREE(svalues);
3901eb62cbbSBarry Smith 
3911eb62cbbSBarry Smith   return 0;
3921eb62cbbSBarry Smith }
3931eb62cbbSBarry Smith 
394416022c9SBarry Smith static int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy)
3951eb62cbbSBarry Smith {
396416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
3971eb62cbbSBarry Smith   int        ierr;
398416022c9SBarry Smith 
39948d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_MPIAIJ:must assemble matrix");
40064119d58SLois Curfman McInnes   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr);
401416022c9SBarry Smith   ierr = MatMult(a->A,xx,yy); CHKERRQ(ierr);
40264119d58SLois Curfman McInnes   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr);
403416022c9SBarry Smith   ierr = MatMultAdd(a->B,a->lvec,yy,yy); CHKERRQ(ierr);
4041eb62cbbSBarry Smith   return 0;
4051eb62cbbSBarry Smith }
4061eb62cbbSBarry Smith 
407416022c9SBarry Smith static int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
408da3a660dSBarry Smith {
409416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
410da3a660dSBarry Smith   int        ierr;
41148d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_MPIAIJ:must assemble matrix");
41264119d58SLois Curfman McInnes   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr);
413416022c9SBarry Smith   ierr = MatMultAdd(a->A,xx,yy,zz); CHKERRQ(ierr);
41464119d58SLois Curfman McInnes   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_ALL,a->Mvctx);CHKERRQ(ierr);
415416022c9SBarry Smith   ierr = MatMultAdd(a->B,a->lvec,zz,zz); CHKERRQ(ierr);
416da3a660dSBarry Smith   return 0;
417da3a660dSBarry Smith }
418da3a660dSBarry Smith 
419416022c9SBarry Smith static int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy)
420da3a660dSBarry Smith {
421416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
422da3a660dSBarry Smith   int        ierr;
423da3a660dSBarry Smith 
42448d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMulTrans_MPIAIJ:must assemble matrix");
425da3a660dSBarry Smith   /* do nondiagonal part */
426416022c9SBarry Smith   ierr = MatMultTrans(a->B,xx,a->lvec); CHKERRQ(ierr);
427da3a660dSBarry Smith   /* send it on its way */
428416022c9SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,
42964119d58SLois Curfman McInnes                 (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr);
430da3a660dSBarry Smith   /* do local part */
431416022c9SBarry Smith   ierr = MatMultTrans(a->A,xx,yy); CHKERRQ(ierr);
432da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
433da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
434da3a660dSBarry Smith   /* but is not perhaps always true. */
435416022c9SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,
43664119d58SLois Curfman McInnes                   (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr);
437da3a660dSBarry Smith   return 0;
438da3a660dSBarry Smith }
439da3a660dSBarry Smith 
440416022c9SBarry Smith static int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
441da3a660dSBarry Smith {
442416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
443da3a660dSBarry Smith   int        ierr;
444da3a660dSBarry Smith 
44548d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMulTransAdd_MPIAIJ:must assemble matrix");
446da3a660dSBarry Smith   /* do nondiagonal part */
447416022c9SBarry Smith   ierr = MatMultTrans(a->B,xx,a->lvec); CHKERRQ(ierr);
448da3a660dSBarry Smith   /* send it on its way */
449416022c9SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,
45064119d58SLois Curfman McInnes                  (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr);
451da3a660dSBarry Smith   /* do local part */
452416022c9SBarry Smith   ierr = MatMultTransAdd(a->A,xx,yy,zz); CHKERRQ(ierr);
453da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
454da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
455da3a660dSBarry Smith   /* but is not perhaps always true. */
456416022c9SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,
45764119d58SLois Curfman McInnes                   (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr);
458da3a660dSBarry Smith   return 0;
459da3a660dSBarry Smith }
460da3a660dSBarry Smith 
4611eb62cbbSBarry Smith /*
4621eb62cbbSBarry Smith   This only works correctly for square matrices where the subblock A->A is the
4631eb62cbbSBarry Smith    diagonal block
4641eb62cbbSBarry Smith */
465416022c9SBarry Smith static int MatGetDiagonal_MPIAIJ(Mat A,Vec v)
4661eb62cbbSBarry Smith {
467416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
46848d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiag_MPIAIJ:must assemble matrix");
469416022c9SBarry Smith   return MatGetDiagonal(a->A,v);
4701eb62cbbSBarry Smith }
4711eb62cbbSBarry Smith 
47244a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj)
4731eb62cbbSBarry Smith {
4741eb62cbbSBarry Smith   Mat        mat = (Mat) obj;
47544a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
4761eb62cbbSBarry Smith   int        ierr;
477a5a9c739SBarry Smith #if defined(PETSC_LOG)
4786e35fa57SLois Curfman McInnes   PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N);
479a5a9c739SBarry Smith #endif
48078b31e54SBarry Smith   PETSCFREE(aij->rowners);
48178b31e54SBarry Smith   ierr = MatDestroy(aij->A); CHKERRQ(ierr);
48278b31e54SBarry Smith   ierr = MatDestroy(aij->B); CHKERRQ(ierr);
48378b31e54SBarry Smith   if (aij->colmap) PETSCFREE(aij->colmap);
48478b31e54SBarry Smith   if (aij->garray) PETSCFREE(aij->garray);
4851eb62cbbSBarry Smith   if (aij->lvec)   VecDestroy(aij->lvec);
4861eb62cbbSBarry Smith   if (aij->Mvctx)  VecScatterCtxDestroy(aij->Mvctx);
48778b31e54SBarry Smith   PETSCFREE(aij);
488a5a9c739SBarry Smith   PLogObjectDestroy(mat);
489a5a9c739SBarry Smith   PETSCHEADERDESTROY(mat);
4901eb62cbbSBarry Smith   return 0;
4911eb62cbbSBarry Smith }
4927857610eSBarry Smith #include "draw.h"
493b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
494ee50ffe9SBarry Smith 
495416022c9SBarry Smith static int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer)
4961eb62cbbSBarry Smith {
497416022c9SBarry Smith   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
498416022c9SBarry Smith   int         ierr;
499416022c9SBarry Smith 
500416022c9SBarry Smith   if (aij->numtids == 1) {
501416022c9SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
502416022c9SBarry Smith   }
503416022c9SBarry Smith   return 0;
504416022c9SBarry Smith }
505416022c9SBarry Smith 
506416022c9SBarry Smith static int MatView_MPIAIJ_ASCIIorDraw(Mat mat,Viewer viewer)
507416022c9SBarry Smith {
50844a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
509dbb450caSBarry Smith   Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data;
510416022c9SBarry Smith   int         ierr, format,shift = C->indexshift;
511d13ab20cSBarry Smith   PetscObject vobj = (PetscObject) viewer;
512d636dbe3SBarry Smith   FILE        *fd;
513416022c9SBarry Smith 
514416022c9SBarry Smith   if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) {
515416022c9SBarry Smith     ierr = ViewerFileGetFormat_Private(viewer,&format);
516416022c9SBarry Smith     if (format == FILE_FORMAT_INFO) {
517416022c9SBarry Smith       return 0;
518416022c9SBarry Smith     }
519416022c9SBarry Smith   }
520416022c9SBarry Smith 
521416022c9SBarry Smith   if (vobj->type == ASCII_FILE_VIEWER) {
522d636dbe3SBarry Smith     ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
5237f813858SBarry Smith     MPIU_Seq_begin(mat->comm,1);
524d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
5251eb62cbbSBarry Smith            aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
5261eb62cbbSBarry Smith            aij->cend);
52778b31e54SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
52878b31e54SBarry Smith     ierr = MatView(aij->B,viewer); CHKERRQ(ierr);
529d13ab20cSBarry Smith     fflush(fd);
5307f813858SBarry Smith     MPIU_Seq_end(mat->comm,1);
531d13ab20cSBarry Smith   }
532416022c9SBarry Smith   else {
53395373324SBarry Smith     int numtids = aij->numtids, mytid = aij->mytid;
53495373324SBarry Smith     if (numtids == 1) {
53578b31e54SBarry Smith       ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
53695373324SBarry Smith     }
53795373324SBarry Smith     else {
53895373324SBarry Smith       /* assemble the entire matrix onto first processor. */
53995373324SBarry Smith       Mat         A;
540ec8511deSBarry Smith       Mat_SeqAIJ *Aloc;
5412eb8c8abSBarry Smith       int         M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
54295373324SBarry Smith       Scalar      *a;
5432ee70a88SLois Curfman McInnes 
54495373324SBarry Smith       if (!mytid) {
545416022c9SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,0,0,0,&A);CHKERRQ(ierr);
54695373324SBarry Smith       }
54795373324SBarry Smith       else {
548416022c9SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,0,0,0,&A);CHKERRQ(ierr);
54995373324SBarry Smith       }
550464493b3SBarry Smith       PLogObjectParent(mat,A);
551416022c9SBarry Smith 
55295373324SBarry Smith 
55395373324SBarry Smith       /* copy over the A part */
554ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->A->data;
5552ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
55695373324SBarry Smith       row = aij->rstart;
557dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;}
55895373324SBarry Smith       for ( i=0; i<m; i++ ) {
559416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr);
56095373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
56195373324SBarry Smith       }
5622ee70a88SLois Curfman McInnes       aj = Aloc->j;
563dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;}
56495373324SBarry Smith 
56595373324SBarry Smith       /* copy over the B part */
566ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->B->data;
5672ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
56895373324SBarry Smith       row = aij->rstart;
56978b31e54SBarry Smith       ct = cols = (int *) PETSCMALLOC( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
570dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];}
57195373324SBarry Smith       for ( i=0; i<m; i++ ) {
572416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr);
57395373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
57495373324SBarry Smith       }
57578b31e54SBarry Smith       PETSCFREE(ct);
57678b31e54SBarry Smith       ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
57778b31e54SBarry Smith       ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
57895373324SBarry Smith       if (!mytid) {
57978b31e54SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr);
58095373324SBarry Smith       }
58178b31e54SBarry Smith       ierr = MatDestroy(A); CHKERRQ(ierr);
58295373324SBarry Smith     }
58395373324SBarry Smith   }
5841eb62cbbSBarry Smith   return 0;
5851eb62cbbSBarry Smith }
5861eb62cbbSBarry Smith 
587416022c9SBarry Smith static int MatView_MPIAIJ(PetscObject obj,Viewer viewer)
588416022c9SBarry Smith {
589416022c9SBarry Smith   Mat         mat = (Mat) obj;
590416022c9SBarry Smith   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
591416022c9SBarry Smith   int         ierr;
592416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
593416022c9SBarry Smith 
59448d91487SBarry Smith   if (!aij->assembled) SETERRQ(1,"MatView_MPIAIJ:must assemble matrix");
595416022c9SBarry Smith   if (!viewer) {
596416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
597416022c9SBarry Smith   }
598416022c9SBarry Smith   if (vobj->cookie == DRAW_COOKIE && vobj->type == NULLWINDOW) return 0;
599416022c9SBarry Smith   else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILE_VIEWER) {
600416022c9SBarry Smith     ierr = MatView_MPIAIJ_ASCIIorDraw(mat,viewer); CHKERRQ(ierr);
601416022c9SBarry Smith   }
602416022c9SBarry Smith   else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILES_VIEWER) {
603416022c9SBarry Smith     ierr = MatView_MPIAIJ_ASCIIorDraw(mat,viewer); CHKERRQ(ierr);
604416022c9SBarry Smith   }
605416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
606416022c9SBarry Smith     ierr = MatView_MPIAIJ_ASCIIorDraw(mat,viewer); CHKERRQ(ierr);
607416022c9SBarry Smith   }
608416022c9SBarry Smith   else if (vobj->type == BINARY_FILE_VIEWER) {
609416022c9SBarry Smith     return MatView_MPIAIJ_Binary(mat,viewer);
610416022c9SBarry Smith   }
611416022c9SBarry Smith   return 0;
612416022c9SBarry Smith }
613416022c9SBarry Smith 
614ec8511deSBarry Smith extern int MatMarkDiag_SeqAIJ(Mat);
6151eb62cbbSBarry Smith /*
6161eb62cbbSBarry Smith     This has to provide several versions.
6171eb62cbbSBarry Smith 
6181eb62cbbSBarry Smith      1) per sequential
6191eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
6201eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
621d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
6221eb62cbbSBarry Smith */
623fc76ce87SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
624dbb450caSBarry Smith                            double fshift,int its,Vec xx)
6258a729477SBarry Smith {
62644a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
627d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
628ec8511deSBarry Smith   Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data;
6296abc6512SBarry Smith   Scalar     zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts;
6306abc6512SBarry Smith   int        ierr,*idx, *diag;
631416022c9SBarry Smith   int        n = mat->n, m = mat->m, i,shift = A->indexshift;
632da3a660dSBarry Smith   Vec        tt;
6338a729477SBarry Smith 
63448d91487SBarry Smith   if (!mat->assembled) SETERRQ(1,"MatRelax_MPIAIJ:must assemble matrix");
6351eb62cbbSBarry Smith 
636d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
637dbb450caSBarry Smith   xs = x + shift; /* shift by one for index start of 1 */
638dbb450caSBarry Smith   ls = ls + shift;
639ec8511deSBarry Smith   if (!A->diag) {if ((ierr = MatMarkDiag_SeqAIJ(AA))) return ierr;}
640d6dfbf8fSBarry Smith   diag = A->diag;
641acb40c82SBarry Smith   if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) {
64248d91487SBarry Smith     SETERRQ(1,"MatRelax_MPIAIJ:Option not supported");
643acb40c82SBarry Smith   }
644da3a660dSBarry Smith   if (flag & SOR_EISENSTAT) {
645da3a660dSBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
646da3a660dSBarry Smith     U is upper triangular, E is diagonal; This routine applies
647da3a660dSBarry Smith 
648da3a660dSBarry Smith             (L + E)^{-1} A (U + E)^{-1}
649da3a660dSBarry Smith 
650da3a660dSBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
651da3a660dSBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
652da3a660dSBarry Smith     is the relaxation factor.
653da3a660dSBarry Smith     */
65478b31e54SBarry Smith     ierr = VecDuplicate(xx,&tt); CHKERRQ(ierr);
655464493b3SBarry Smith     PLogObjectParent(matin,tt);
656da3a660dSBarry Smith     VecGetArray(tt,&t);
657da3a660dSBarry Smith     scale = (2.0/omega) - 1.0;
658da3a660dSBarry Smith     /*  x = (E + U)^{-1} b */
659da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
66064119d58SLois Curfman McInnes     ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
66178b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
662da3a660dSBarry Smith     for ( i=m-1; i>-1; i-- ) {
663da3a660dSBarry Smith       n    = A->i[i+1] - diag[i] - 1;
664dbb450caSBarry Smith       idx  = A->j + diag[i] + !shift;
665dbb450caSBarry Smith       v    = A->a + diag[i] + !shift;
666da3a660dSBarry Smith       sum  = b[i];
667da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
668dbb450caSBarry Smith       d    = fshift + A->a[diag[i]+shift];
669da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
670dbb450caSBarry Smith       idx  = B->j + B->i[i] + shift;
671dbb450caSBarry Smith       v    = B->a + B->i[i] + shift;
672da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
673da3a660dSBarry Smith       x[i] = omega*(sum/d);
674da3a660dSBarry Smith     }
67564119d58SLois Curfman McInnes     ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
67678b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
677da3a660dSBarry Smith 
678da3a660dSBarry Smith     /*  t = b - (2*E - D)x */
679da3a660dSBarry Smith     v = A->a;
680dbb450caSBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
681da3a660dSBarry Smith 
682da3a660dSBarry Smith     /*  t = (E + L)^{-1}t */
683dbb450caSBarry Smith     ts = t + shift; /* shifted by one for index start of a or mat->j*/
684da3a660dSBarry Smith     diag = A->diag;
685da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
68664119d58SLois Curfman McInnes     ierr = VecPipelineBegin(tt,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
68778b31e54SBarry Smith                                                  mat->Mvctx); CHKERRQ(ierr);
688da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
689da3a660dSBarry Smith       n    = diag[i] - A->i[i];
690dbb450caSBarry Smith       idx  = A->j + A->i[i] + shift;
691dbb450caSBarry Smith       v    = A->a + A->i[i] + shift;
692da3a660dSBarry Smith       sum  = t[i];
693da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
694dbb450caSBarry Smith       d    = fshift + A->a[diag[i]+shift];
695da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
696dbb450caSBarry Smith       idx  = B->j + B->i[i] + shift;
697dbb450caSBarry Smith       v    = B->a + B->i[i] + shift;
698da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
699da3a660dSBarry Smith       t[i] = omega*(sum/d);
700da3a660dSBarry Smith     }
70164119d58SLois Curfman McInnes     ierr = VecPipelineEnd(tt,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
70278b31e54SBarry Smith                                                     mat->Mvctx); CHKERRQ(ierr);
703da3a660dSBarry Smith     /*  x = x + t */
704da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
705da3a660dSBarry Smith     VecDestroy(tt);
706da3a660dSBarry Smith     return 0;
707da3a660dSBarry Smith   }
708da3a660dSBarry Smith 
7091eb62cbbSBarry Smith 
710d6dfbf8fSBarry Smith   if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){
711da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
712da3a660dSBarry Smith       VecSet(&zero,mat->lvec); VecSet(&zero,xx);
713da3a660dSBarry Smith     }
714da3a660dSBarry Smith     else {
71564119d58SLois Curfman McInnes       ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx);
71678b31e54SBarry Smith       CHKERRQ(ierr);
71764119d58SLois Curfman McInnes       ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx);
71878b31e54SBarry Smith       CHKERRQ(ierr);
719da3a660dSBarry Smith     }
720d6dfbf8fSBarry Smith     while (its--) {
721d6dfbf8fSBarry Smith       /* go down through the rows */
72264119d58SLois Curfman McInnes       ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
72378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
724d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
725d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
726dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
727dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
728d6dfbf8fSBarry Smith         sum  = b[i];
729d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
730dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
731d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
732dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
733dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
734d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
735d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
736d6dfbf8fSBarry Smith       }
73764119d58SLois Curfman McInnes       ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
73878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
739d6dfbf8fSBarry Smith       /* come up through the rows */
74064119d58SLois Curfman McInnes       ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
74178b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
742d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
743d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
744dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
745dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
746d6dfbf8fSBarry Smith         sum  = b[i];
747d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
748dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
749d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
750dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
751dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
752d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
753d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
754d6dfbf8fSBarry Smith       }
75564119d58SLois Curfman McInnes       ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
75678b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
757d6dfbf8fSBarry Smith     }
758d6dfbf8fSBarry Smith   }
759d6dfbf8fSBarry Smith   else if (flag & SOR_FORWARD_SWEEP){
760da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
761da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
76264119d58SLois Curfman McInnes       ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
76378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
764da3a660dSBarry Smith       for ( i=0; i<m; i++ ) {
765da3a660dSBarry Smith         n    = diag[i] - A->i[i];
766dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
767dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
768da3a660dSBarry Smith         sum  = b[i];
769da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
770dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
771da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
772dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
773dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
774da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
775da3a660dSBarry Smith         x[i] = omega*(sum/d);
776da3a660dSBarry Smith       }
77764119d58SLois Curfman McInnes       ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
77878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
779da3a660dSBarry Smith       its--;
780da3a660dSBarry Smith     }
781d6dfbf8fSBarry Smith     while (its--) {
78264119d58SLois Curfman McInnes       ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx);
78378b31e54SBarry Smith       CHKERRQ(ierr);
78464119d58SLois Curfman McInnes       ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_UP,mat->Mvctx);
78578b31e54SBarry Smith       CHKERRQ(ierr);
78664119d58SLois Curfman McInnes       ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
78778b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
788d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
789d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
790dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
791dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
792d6dfbf8fSBarry Smith         sum  = b[i];
793d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
794dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
795d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
796dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
797dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
798d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
799d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
800d6dfbf8fSBarry Smith       }
80164119d58SLois Curfman McInnes       ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_DOWN,
80278b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
803d6dfbf8fSBarry Smith     }
804d6dfbf8fSBarry Smith   }
805d6dfbf8fSBarry Smith   else if (flag & SOR_BACKWARD_SWEEP){
806da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
807da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
80864119d58SLois Curfman McInnes       ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
80978b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
810da3a660dSBarry Smith       for ( i=m-1; i>-1; i-- ) {
811da3a660dSBarry Smith         n    = A->i[i+1] - diag[i] - 1;
812dbb450caSBarry Smith         idx  = A->j + diag[i] + !shift;
813dbb450caSBarry Smith         v    = A->a + diag[i] + !shift;
814da3a660dSBarry Smith         sum  = b[i];
815da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
816dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
817da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
818dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
819dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
820da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
821da3a660dSBarry Smith         x[i] = omega*(sum/d);
822da3a660dSBarry Smith       }
82364119d58SLois Curfman McInnes       ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
82478b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
825da3a660dSBarry Smith       its--;
826da3a660dSBarry Smith     }
827d6dfbf8fSBarry Smith     while (its--) {
82864119d58SLois Curfman McInnes       ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_DOWN,
82978b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
83064119d58SLois Curfman McInnes       ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_DOWN,
83178b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
83264119d58SLois Curfman McInnes       ierr = VecPipelineBegin(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
83378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
834d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
835d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
836dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
837dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
838d6dfbf8fSBarry Smith         sum  = b[i];
839d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
840dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
841d6dfbf8fSBarry 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;
844d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
845d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
846d6dfbf8fSBarry Smith       }
84764119d58SLois Curfman McInnes       ierr = VecPipelineEnd(xx,mat->lvec,INSERT_VALUES,PIPELINE_UP,
84878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
849d6dfbf8fSBarry Smith     }
850d6dfbf8fSBarry Smith   }
851d6dfbf8fSBarry Smith   else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
852da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
853dbb450caSBarry Smith       return MatRelax(mat->A,bb,omega,flag,fshift,its,xx);
854da3a660dSBarry Smith     }
85564119d58SLois Curfman McInnes     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx);
85678b31e54SBarry Smith     CHKERRQ(ierr);
85764119d58SLois Curfman McInnes     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx);
85878b31e54SBarry Smith     CHKERRQ(ierr);
859d6dfbf8fSBarry Smith     while (its--) {
860d6dfbf8fSBarry Smith       /* go down through the rows */
861d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
862d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
863dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
864dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
865d6dfbf8fSBarry Smith         sum  = b[i];
866d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
867dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
868d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
869dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
870dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
871d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
872d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
873d6dfbf8fSBarry Smith       }
874d6dfbf8fSBarry Smith       /* come up through the rows */
875d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
876d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
877dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
878dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
879d6dfbf8fSBarry Smith         sum  = b[i];
880d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
881dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
882d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
883dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
884dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
885d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
886d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
887d6dfbf8fSBarry Smith       }
888d6dfbf8fSBarry Smith     }
889d6dfbf8fSBarry Smith   }
890d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_FORWARD_SWEEP){
891da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
892dbb450caSBarry Smith       return MatRelax(mat->A,bb,omega,flag,fshift,its,xx);
893da3a660dSBarry Smith     }
89464119d58SLois Curfman McInnes     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx);
89578b31e54SBarry Smith     CHKERRQ(ierr);
89664119d58SLois Curfman McInnes     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,mat->Mvctx);
89778b31e54SBarry Smith     CHKERRQ(ierr);
898d6dfbf8fSBarry Smith     while (its--) {
899d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
900d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
901dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
902dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
903d6dfbf8fSBarry Smith         sum  = b[i];
904d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
905dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
906d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
907dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
908dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
909d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
910d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
911d6dfbf8fSBarry Smith       }
912d6dfbf8fSBarry Smith     }
913d6dfbf8fSBarry Smith   }
914d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_BACKWARD_SWEEP){
915da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
916dbb450caSBarry Smith       return MatRelax(mat->A,bb,omega,flag,fshift,its,xx);
917da3a660dSBarry Smith     }
91864119d58SLois Curfman McInnes     ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,
91978b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
92064119d58SLois Curfman McInnes     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_ALL,
92178b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
922d6dfbf8fSBarry Smith     while (its--) {
923d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
924d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
925dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
926dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
927d6dfbf8fSBarry Smith         sum  = b[i];
928d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
929dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
930d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
931dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
932dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
933d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
934d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
935d6dfbf8fSBarry Smith       }
936d6dfbf8fSBarry Smith     }
937d6dfbf8fSBarry Smith   }
9388a729477SBarry Smith   return 0;
9398a729477SBarry Smith }
940a66be287SLois Curfman McInnes 
941fc76ce87SLois Curfman McInnes static int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,int *nz,
942fc76ce87SLois Curfman McInnes                              int *nzalloc,int *mem)
943a66be287SLois Curfman McInnes {
944a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
945a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
946a66be287SLois Curfman McInnes   int        ierr, isend[3], irecv[3], nzA, nzallocA, memA;
947a66be287SLois Curfman McInnes 
94878b31e54SBarry Smith   ierr = MatGetInfo(A,MAT_LOCAL,&nzA,&nzallocA,&memA); CHKERRQ(ierr);
94978b31e54SBarry Smith   ierr = MatGetInfo(B,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr);
950a66be287SLois Curfman McInnes   isend[0] += nzA; isend[1] += nzallocA; isend[2] += memA;
951a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
952a66be287SLois Curfman McInnes     *nz = isend[0]; *nzalloc = isend[1]; *mem = isend[2];
953a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
954d65a2f8fSBarry Smith     MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_MAX,matin->comm);
955a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
956a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
957d65a2f8fSBarry Smith     MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_SUM,matin->comm);
958a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
959a66be287SLois Curfman McInnes   }
960a66be287SLois Curfman McInnes   return 0;
961a66be287SLois Curfman McInnes }
962a66be287SLois Curfman McInnes 
963299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*);
964299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*);
965299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double);
966299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *);
967299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec);
968299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec);
969299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec);
970299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec);
971299609e3SLois Curfman McInnes 
972416022c9SBarry Smith static int MatSetOption_MPIAIJ(Mat A,MatOption op)
973c74985f6SBarry Smith {
974c0bbcb79SLois Curfman McInnes   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
975c74985f6SBarry Smith 
976c0bbcb79SLois Curfman McInnes   if (op == NO_NEW_NONZERO_LOCATIONS ||
977c0bbcb79SLois Curfman McInnes       op == YES_NEW_NONZERO_LOCATIONS ||
978c0bbcb79SLois Curfman McInnes       op == COLUMNS_SORTED ||
979c0bbcb79SLois Curfman McInnes       op == ROW_ORIENTED) {
980c0bbcb79SLois Curfman McInnes         MatSetOption(a->A,op);
981c0bbcb79SLois Curfman McInnes         MatSetOption(a->B,op);
982c74985f6SBarry Smith   }
983c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
984c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
985c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
986c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
987c0bbcb79SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_MPIAIJ:Option ignored\n");
988bbb6d6a8SBarry Smith   else if (op == COLUMN_ORIENTED)
989*4d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:COLUMN_ORIENTED");}
990c0bbcb79SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
991*4d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:NO_NEW_DIAGONALS");}
992c0bbcb79SLois Curfman McInnes   else
993*4d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIAIJ:unknown option");}
994c74985f6SBarry Smith   return 0;
995c74985f6SBarry Smith }
996c74985f6SBarry Smith 
997d1710a03SLois Curfman McInnes static int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
998c74985f6SBarry Smith {
99944a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1000c74985f6SBarry Smith   *m = mat->M; *n = mat->N;
1001c74985f6SBarry Smith   return 0;
1002c74985f6SBarry Smith }
1003c74985f6SBarry Smith 
1004d1710a03SLois Curfman McInnes static int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
1005c74985f6SBarry Smith {
100644a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1007b7c46309SBarry Smith   *m = mat->m; *n = mat->N;
1008c74985f6SBarry Smith   return 0;
1009c74985f6SBarry Smith }
1010c74985f6SBarry Smith 
1011d1710a03SLois Curfman McInnes static int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
1012c74985f6SBarry Smith {
101344a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1014c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
1015c74985f6SBarry Smith   return 0;
1016c74985f6SBarry Smith }
1017c74985f6SBarry Smith 
101839e00950SLois Curfman McInnes static int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
101939e00950SLois Curfman McInnes {
1020154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1021154123eaSLois Curfman McInnes   Scalar     *vworkA, *vworkB, **pvA, **pvB;
1022154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
1023154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
102439e00950SLois Curfman McInnes 
1025416022c9SBarry Smith   if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIAIJ:only local rows")
1026abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
102739e00950SLois Curfman McInnes 
1028154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1029154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1030154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
103178b31e54SBarry Smith   ierr = MatGetRow(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
103278b31e54SBarry Smith   ierr = MatGetRow(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
1033154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1034154123eaSLois Curfman McInnes 
1035154123eaSLois Curfman McInnes   if (v  || idx) {
1036154123eaSLois Curfman McInnes     if (nztot) {
1037154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
1038299609e3SLois Curfman McInnes       int imark;
103948b35521SBarry Smith       if (mat->assembled) {
1040154123eaSLois Curfman McInnes         for (i=0; i<nzB; i++) cworkB[i] = mat->garray[cworkB[i]];
104148b35521SBarry Smith       }
1042154123eaSLois Curfman McInnes       if (v) {
104378b31e54SBarry Smith         *v = (Scalar *) PETSCMALLOC( (nztot)*sizeof(Scalar) ); CHKPTRQ(*v);
104439e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1045154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*v)[i] = vworkB[i];
1046154123eaSLois Curfman McInnes           else break;
1047154123eaSLois Curfman McInnes         }
1048154123eaSLois Curfman McInnes         imark = i;
1049154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*v)[imark+i] = vworkA[i];
1050299609e3SLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*v)[nzA+i] = vworkB[i];
1051154123eaSLois Curfman McInnes       }
1052154123eaSLois Curfman McInnes       if (idx) {
105378b31e54SBarry Smith         *idx = (int *) PETSCMALLOC( (nztot)*sizeof(int) ); CHKPTRQ(*idx);
1054154123eaSLois Curfman McInnes         for (i=0; i<nzA; i++) cworkA[i] += cstart;
1055154123eaSLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1056154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*idx)[i] = cworkB[i];
1057154123eaSLois Curfman McInnes           else break;
1058154123eaSLois Curfman McInnes         }
1059154123eaSLois Curfman McInnes         imark = i;
1060154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*idx)[imark+i] = cworkA[i];
1061299609e3SLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*idx)[nzA+i] = cworkB[i];
106239e00950SLois Curfman McInnes       }
106339e00950SLois Curfman McInnes     }
106439e00950SLois Curfman McInnes     else {*idx = 0; *v=0;}
1065154123eaSLois Curfman McInnes   }
106639e00950SLois Curfman McInnes   *nz = nztot;
106778b31e54SBarry Smith   ierr = MatRestoreRow(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
106878b31e54SBarry Smith   ierr = MatRestoreRow(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
106939e00950SLois Curfman McInnes   return 0;
107039e00950SLois Curfman McInnes }
107139e00950SLois Curfman McInnes 
107239e00950SLois Curfman McInnes static int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
107339e00950SLois Curfman McInnes {
107478b31e54SBarry Smith   if (idx) PETSCFREE(*idx);
107578b31e54SBarry Smith   if (v) PETSCFREE(*v);
107639e00950SLois Curfman McInnes   return 0;
107739e00950SLois Curfman McInnes }
107839e00950SLois Curfman McInnes 
1079855ac2c5SLois Curfman McInnes static int MatNorm_MPIAIJ(Mat mat,MatNormType type,double *norm)
1080855ac2c5SLois Curfman McInnes {
1081855ac2c5SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1082ec8511deSBarry Smith   Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data;
1083416022c9SBarry Smith   int        ierr, i, j, cstart = aij->cstart,shift = amat->indexshift;
1084416022c9SBarry Smith   double     sum = 0.0;
108504ca555eSLois Curfman McInnes   Scalar     *v;
108604ca555eSLois Curfman McInnes 
1087416022c9SBarry Smith   if (!aij->assembled) SETERRQ(1,"MatNorm_MPIAIJ:Must assemble matrix");
1088855ac2c5SLois Curfman McInnes   if (aij->numtids == 1) {
108914183eadSLois Curfman McInnes     ierr =  MatNorm(aij->A,type,norm); CHKERRQ(ierr);
109037fa93a5SLois Curfman McInnes   } else {
109104ca555eSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
109204ca555eSLois Curfman McInnes       v = amat->a;
109304ca555eSLois Curfman McInnes       for (i=0; i<amat->nz; i++ ) {
109404ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
109504ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
109604ca555eSLois Curfman McInnes #else
109704ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
109804ca555eSLois Curfman McInnes #endif
109904ca555eSLois Curfman McInnes       }
110004ca555eSLois Curfman McInnes       v = bmat->a;
110104ca555eSLois Curfman McInnes       for (i=0; i<bmat->nz; i++ ) {
110204ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
110304ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
110404ca555eSLois Curfman McInnes #else
110504ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
110604ca555eSLois Curfman McInnes #endif
110704ca555eSLois Curfman McInnes       }
110804ca555eSLois Curfman McInnes       MPI_Allreduce((void*)&sum,(void*)norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);
110904ca555eSLois Curfman McInnes       *norm = sqrt(*norm);
111004ca555eSLois Curfman McInnes     }
111104ca555eSLois Curfman McInnes     else if (type == NORM_1) { /* max column norm */
111204ca555eSLois Curfman McInnes       double *tmp, *tmp2;
111304ca555eSLois Curfman McInnes       int    *jj, *garray = aij->garray;
111404ca555eSLois Curfman McInnes       tmp  = (double *) PETSCMALLOC( aij->N*sizeof(double) ); CHKPTRQ(tmp);
111504ca555eSLois Curfman McInnes       tmp2 = (double *) PETSCMALLOC( aij->N*sizeof(double) ); CHKPTRQ(tmp2);
1116416022c9SBarry Smith       PetscZero(tmp,aij->N*sizeof(double));
111704ca555eSLois Curfman McInnes       *norm = 0.0;
111804ca555eSLois Curfman McInnes       v = amat->a; jj = amat->j;
111904ca555eSLois Curfman McInnes       for ( j=0; j<amat->nz; j++ ) {
112004ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
1121dbb450caSBarry Smith         tmp[cstart + *jj++ + shift] += abs(*v++);
112204ca555eSLois Curfman McInnes #else
1123dbb450caSBarry Smith         tmp[cstart + *jj++ + shift] += fabs(*v++);
112404ca555eSLois Curfman McInnes #endif
112504ca555eSLois Curfman McInnes       }
112604ca555eSLois Curfman McInnes       v = bmat->a; jj = bmat->j;
112704ca555eSLois Curfman McInnes       for ( j=0; j<bmat->nz; j++ ) {
112804ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
1129dbb450caSBarry Smith         tmp[garray[*jj++ + shift]] += abs(*v++);
113004ca555eSLois Curfman McInnes #else
1131dbb450caSBarry Smith         tmp[garray[*jj++ + shift]] += fabs(*v++);
113204ca555eSLois Curfman McInnes #endif
113304ca555eSLois Curfman McInnes       }
113404ca555eSLois Curfman McInnes       MPI_Allreduce((void*)tmp,(void*)tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);
113504ca555eSLois Curfman McInnes       for ( j=0; j<aij->N; j++ ) {
113604ca555eSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
113704ca555eSLois Curfman McInnes       }
113804ca555eSLois Curfman McInnes       PETSCFREE(tmp); PETSCFREE(tmp2);
113904ca555eSLois Curfman McInnes     }
114004ca555eSLois Curfman McInnes     else if (type == NORM_INFINITY) { /* max row norm */
114104ca555eSLois Curfman McInnes       double normtemp = 0.0;
114204ca555eSLois Curfman McInnes       for ( j=0; j<amat->m; j++ ) {
1143dbb450caSBarry Smith         v = amat->a + amat->i[j] + shift;
114404ca555eSLois Curfman McInnes         sum = 0.0;
114504ca555eSLois Curfman McInnes         for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) {
114604ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
114704ca555eSLois Curfman McInnes           sum += abs(*v); v++;
114804ca555eSLois Curfman McInnes #else
114904ca555eSLois Curfman McInnes           sum += fabs(*v); v++;
115004ca555eSLois Curfman McInnes #endif
115104ca555eSLois Curfman McInnes         }
1152dbb450caSBarry Smith         v = bmat->a + bmat->i[j] + shift;
115304ca555eSLois Curfman McInnes         for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) {
115404ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
115504ca555eSLois Curfman McInnes           sum += abs(*v); v++;
115604ca555eSLois Curfman McInnes #else
115704ca555eSLois Curfman McInnes           sum += fabs(*v); v++;
115804ca555eSLois Curfman McInnes #endif
115904ca555eSLois Curfman McInnes         }
116004ca555eSLois Curfman McInnes         if (sum > normtemp) normtemp = sum;
116104ca555eSLois Curfman McInnes         MPI_Allreduce((void*)&normtemp,(void*)norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);
116204ca555eSLois Curfman McInnes       }
116304ca555eSLois Curfman McInnes     }
116404ca555eSLois Curfman McInnes     else {
116548d91487SBarry Smith       SETERRQ(1,"MatNorm_MPIRow:No support for two norm");
116604ca555eSLois Curfman McInnes     }
116737fa93a5SLois Curfman McInnes   }
1168855ac2c5SLois Curfman McInnes   return 0;
1169855ac2c5SLois Curfman McInnes }
1170855ac2c5SLois Curfman McInnes 
11710de55854SLois Curfman McInnes static int MatTranspose_MPIAIJ(Mat A,Mat *matout)
1172b7c46309SBarry Smith {
1173b7c46309SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1174dbb450caSBarry Smith   Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data;
1175416022c9SBarry Smith   int        ierr,shift = Aloc->indexshift;
1176416022c9SBarry Smith   Mat        B;
1177b7c46309SBarry Smith   int        M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct;
1178b7c46309SBarry Smith   Scalar     *array;
1179b7c46309SBarry Smith 
1180416022c9SBarry Smith   if (!matout && M != N) SETERRQ(1,"MatTranspose_MPIAIJ:Square only in-place");
1181b7c46309SBarry Smith   ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,0,0,0,&B);
1182b7c46309SBarry Smith   CHKERRQ(ierr);
1183b7c46309SBarry Smith 
1184b7c46309SBarry Smith   /* copy over the A part */
1185ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->A->data;
1186b7c46309SBarry Smith   m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1187b7c46309SBarry Smith   row = a->rstart;
1188dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;}
1189b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1190416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1191b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
1192b7c46309SBarry Smith   }
1193b7c46309SBarry Smith   aj = Aloc->j;
1194dbb450caSBarry Smith   for ( i=0; i<ai[m]|+shift; i++ ) {aj[i] -= a->cstart + shift;}
1195b7c46309SBarry Smith 
1196b7c46309SBarry Smith   /* copy over the B part */
1197ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->B->data;
1198b7c46309SBarry Smith   m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1199b7c46309SBarry Smith   row = a->rstart;
12001987afe7SBarry Smith   ct = cols = (int *) PETSCMALLOC( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols);
1201dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];}
1202b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1203416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1204b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
1205b7c46309SBarry Smith   }
1206b7c46309SBarry Smith   PETSCFREE(ct);
1207b7c46309SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1208b7c46309SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
12090de55854SLois Curfman McInnes   if (matout) {
12100de55854SLois Curfman McInnes     *matout = B;
12110de55854SLois Curfman McInnes   } else {
12120de55854SLois Curfman McInnes     /* This isn't really an in-place transpose .... but free data structures from a */
12130de55854SLois Curfman McInnes     PETSCFREE(a->rowners);
12140de55854SLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
12150de55854SLois Curfman McInnes     ierr = MatDestroy(a->B); CHKERRQ(ierr);
12160de55854SLois Curfman McInnes     if (a->colmap) PETSCFREE(a->colmap);
12170de55854SLois Curfman McInnes     if (a->garray) PETSCFREE(a->garray);
12180de55854SLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
12190de55854SLois Curfman McInnes     if (a->Mvctx) VecScatterCtxDestroy(a->Mvctx);
12200de55854SLois Curfman McInnes     PETSCFREE(a);
1221416022c9SBarry Smith     PetscMemcpy(A,B,sizeof(struct _Mat));
12220de55854SLois Curfman McInnes     PETSCHEADERDESTROY(B);
12230de55854SLois Curfman McInnes   }
1224b7c46309SBarry Smith   return 0;
1225b7c46309SBarry Smith }
1226b7c46309SBarry Smith 
1227fc76ce87SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MatType,Mat *);
1228ff7509f2SBarry Smith static int MatCopyPrivate_MPIAIJ(Mat,Mat *);
1229d6dfbf8fSBarry Smith 
12308a729477SBarry Smith /* -------------------------------------------------------------------*/
12312ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
123239e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
123344a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
123444a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
1235299609e3SLois Curfman McInnes        MatSolve_MPIAIJ,MatSolveAdd_MPIAIJ,
1236299609e3SLois Curfman McInnes        MatSolveTrans_MPIAIJ,MatSolveTransAdd_MPIAIJ,
1237299609e3SLois Curfman McInnes        MatLUFactor_MPIAIJ,0,
123844a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
1239b7c46309SBarry Smith        MatTranspose_MPIAIJ,
1240a66be287SLois Curfman McInnes        MatGetInfo_MPIAIJ,0,
1241855ac2c5SLois Curfman McInnes        MatGetDiagonal_MPIAIJ,0,MatNorm_MPIAIJ,
1242ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
12431eb62cbbSBarry Smith        0,
1244299609e3SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,
1245299609e3SLois Curfman McInnes        MatGetReordering_MPIAIJ,
1246299609e3SLois Curfman McInnes        MatLUFactorSymbolic_MPIAIJ,MatLUFactorNumeric_MPIAIJ,0,0,
1247d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
1248299609e3SLois Curfman McInnes        MatILUFactorSymbolic_MPIAIJ,0,
1249ff7509f2SBarry Smith        0,0,MatConvert_MPIAIJ,0,0,MatCopyPrivate_MPIAIJ};
12508a729477SBarry Smith 
12511987afe7SBarry Smith /*@C
1252ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
1253ff756334SLois Curfman McInnes    (the default parallel PETSc format).
12548a729477SBarry Smith 
12558a729477SBarry Smith    Input Parameters:
12561eb62cbbSBarry Smith .  comm - MPI communicator
12577d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
12587d3e4905SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated
12597d3e4905SLois Curfman McInnes            if N is given)
12607d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
12617d3e4905SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated
12627d3e4905SLois Curfman McInnes            if n is given)
1263ab693e5aSLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of local submatrix
1264ff756334SLois Curfman McInnes            (same for all local rows)
1265ab693e5aSLois Curfman McInnes .  d_nzz - number of nonzeros per row in diagonal portion of local submatrix
1266ab693e5aSLois Curfman McInnes            or null (possibly different for each row).  You must leave room
1267ab693e5aSLois Curfman McInnes            for the diagonal entry even if it is zero.
1268ab693e5aSLois Curfman McInnes .  o_nz - number of nonzeros per row in off-diagonal portion of local
1269ab693e5aSLois Curfman McInnes            submatrix (same for all local rows).
1270ab693e5aSLois Curfman McInnes .  o_nzz - number of nonzeros per row in off-diagonal portion of local
1271ab693e5aSLois Curfman McInnes            submatrix or null (possibly different for each row).
12728a729477SBarry Smith 
1273ff756334SLois Curfman McInnes    Output Parameter:
12748a729477SBarry Smith .  newmat - the matrix
12758a729477SBarry Smith 
1276ff756334SLois Curfman McInnes    Notes:
1277ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1278ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
1279ff756334SLois Curfman McInnes    storage.  That is, the stored row and column indices begin at
1280ab693e5aSLois Curfman McInnes    one, not zero.  See the users manual for further details.
1281ff756334SLois Curfman McInnes 
1282ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1283ff756334SLois Curfman McInnes    (possibly both).
1284ff756334SLois Curfman McInnes 
1285e0245417SLois Curfman McInnes    Storage Information:
1286e0245417SLois Curfman McInnes    For a square global matrix we define each processor's diagonal portion
1287e0245417SLois Curfman McInnes    to be its local rows and the corresponding columns (a square submatrix);
1288e0245417SLois Curfman McInnes    each processor's off-diagonal portion encompasses the remainder of the
1289e0245417SLois Curfman McInnes    local matrix (a rectangular submatrix).
1290e0245417SLois Curfman McInnes 
1291e0245417SLois Curfman McInnes    The user can specify preallocated storage for the diagonal part of
1292e0245417SLois Curfman McInnes    the local submatrix with either d_nz or d_nnz (not both). Set both
1293e0245417SLois Curfman McInnes    d_nz and d_nnz to zero for PETSc to control dynamic memory allocation.
1294e0245417SLois Curfman McInnes    Likewise, specify preallocated storage for the off-diagonal part of
1295e0245417SLois Curfman McInnes    the local submatrix with o_nz or o_nnz (not both).
1296ff756334SLois Curfman McInnes 
1297dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel
1298ff756334SLois Curfman McInnes 
1299fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
13008a729477SBarry Smith @*/
13011eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,
13021eb62cbbSBarry Smith                     int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat)
13038a729477SBarry Smith {
13048a729477SBarry Smith   Mat          mat;
1305416022c9SBarry Smith   Mat_MPIAIJ   *a;
13066abc6512SBarry Smith   int          ierr, i,sum[2],work[2];
1307416022c9SBarry Smith 
13088a729477SBarry Smith   *newmat         = 0;
130944a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm);
1310a5a9c739SBarry Smith   PLogObjectCreate(mat);
1311416022c9SBarry Smith   mat->data       = (void *) (a = PETSCNEW(Mat_MPIAIJ)); CHKPTRQ(a);
1312416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
131344a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
131444a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
13158a729477SBarry Smith   mat->factor     = 0;
1316d6dfbf8fSBarry Smith 
131764119d58SLois Curfman McInnes   a->insertmode = NOT_SET_VALUES;
1318416022c9SBarry Smith   MPI_Comm_rank(comm,&a->mytid);
1319416022c9SBarry Smith   MPI_Comm_size(comm,&a->numtids);
13201eb62cbbSBarry Smith 
13211987afe7SBarry Smith   if (m == PETSC_DECIDE && (d_nnz || o_nnz))
13221987afe7SBarry Smith     SETERRQ(1,"MatCreateMPIAIJ:Cannot have PETSc decide rows but set d_nnz or o_nnz");
13231987afe7SBarry Smith 
1324dbd7a890SLois Curfman McInnes   if (M == PETSC_DECIDE || N == PETSC_DECIDE) {
13251eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1326d65a2f8fSBarry Smith     MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm );
1327dbd7a890SLois Curfman McInnes     if (M == PETSC_DECIDE) M = sum[0];
1328dbd7a890SLois Curfman McInnes     if (N == PETSC_DECIDE) N = sum[1];
13291eb62cbbSBarry Smith   }
1330416022c9SBarry Smith   if (m == PETSC_DECIDE) {m = M/a->numtids + ((M % a->numtids) > a->mytid);}
1331416022c9SBarry Smith   if (n == PETSC_DECIDE) {n = N/a->numtids + ((N % a->numtids) > a->mytid);}
1332416022c9SBarry Smith   a->m = m;
1333416022c9SBarry Smith   a->n = n;
1334416022c9SBarry Smith   a->N = N;
1335416022c9SBarry Smith   a->M = M;
13361eb62cbbSBarry Smith 
13371eb62cbbSBarry Smith   /* build local table of row and column ownerships */
1338416022c9SBarry Smith   a->rowners = (int *) PETSCMALLOC(2*(a->numtids+2)*sizeof(int)); CHKPTRQ(a->rowners);
1339416022c9SBarry Smith   PLogObjectMemory(mat,2*(a->numtids+2)*sizeof(int)+sizeof(struct _Mat)+
1340464493b3SBarry Smith                        sizeof(Mat_MPIAIJ));
1341416022c9SBarry Smith   a->cowners = a->rowners + a->numtids + 1;
1342416022c9SBarry Smith   MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);
1343416022c9SBarry Smith   a->rowners[0] = 0;
1344416022c9SBarry Smith   for ( i=2; i<=a->numtids; i++ ) {
1345416022c9SBarry Smith     a->rowners[i] += a->rowners[i-1];
13468a729477SBarry Smith   }
1347416022c9SBarry Smith   a->rstart = a->rowners[a->mytid];
1348416022c9SBarry Smith   a->rend   = a->rowners[a->mytid+1];
1349416022c9SBarry Smith   MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);
1350416022c9SBarry Smith   a->cowners[0] = 0;
1351416022c9SBarry Smith   for ( i=2; i<=a->numtids; i++ ) {
1352416022c9SBarry Smith     a->cowners[i] += a->cowners[i-1];
13538a729477SBarry Smith   }
1354416022c9SBarry Smith   a->cstart = a->cowners[a->mytid];
1355416022c9SBarry Smith   a->cend   = a->cowners[a->mytid+1];
13568a729477SBarry Smith 
1357416022c9SBarry Smith   ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&a->A); CHKERRQ(ierr);
1358416022c9SBarry Smith   PLogObjectParent(mat,a->A);
1359416022c9SBarry Smith   ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&a->B); CHKERRQ(ierr);
1360416022c9SBarry Smith   PLogObjectParent(mat,a->B);
13618a729477SBarry Smith 
13621eb62cbbSBarry Smith   /* build cache for off array entries formed */
1363416022c9SBarry Smith   ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr);
1364416022c9SBarry Smith   a->colmap    = 0;
1365416022c9SBarry Smith   a->garray    = 0;
13668a729477SBarry Smith 
13671eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
1368416022c9SBarry Smith   a->lvec      = 0;
1369416022c9SBarry Smith   a->Mvctx     = 0;
1370416022c9SBarry Smith   a->assembled = 0;
13718a729477SBarry Smith 
1372d6dfbf8fSBarry Smith   *newmat = mat;
1373d6dfbf8fSBarry Smith   return 0;
1374d6dfbf8fSBarry Smith }
1375c74985f6SBarry Smith 
1376ff7509f2SBarry Smith static int MatCopyPrivate_MPIAIJ(Mat matin,Mat *newmat)
1377d6dfbf8fSBarry Smith {
1378d6dfbf8fSBarry Smith   Mat        mat;
1379416022c9SBarry Smith   Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data;
13802ee70a88SLois Curfman McInnes   int        ierr, len;
1381d6dfbf8fSBarry Smith 
1382416022c9SBarry Smith   if (!oldmat->assembled) SETERRQ(1,"MatCopyPrivate_MPIAIJ:Must assemble matrix");
1383416022c9SBarry Smith   *newmat      = 0;
138444a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm);
1385a5a9c739SBarry Smith   PLogObjectCreate(mat);
1386416022c9SBarry Smith   mat->data       = (void *) (a = PETSCNEW(Mat_MPIAIJ)); CHKPTRQ(a);
1387416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
138844a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
138944a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
1390d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1391d6dfbf8fSBarry Smith 
1392416022c9SBarry Smith   a->m          = oldmat->m;
1393416022c9SBarry Smith   a->n          = oldmat->n;
1394416022c9SBarry Smith   a->M          = oldmat->M;
1395416022c9SBarry Smith   a->N          = oldmat->N;
1396d6dfbf8fSBarry Smith 
1397416022c9SBarry Smith   a->assembled  = 1;
1398416022c9SBarry Smith   a->rstart     = oldmat->rstart;
1399416022c9SBarry Smith   a->rend       = oldmat->rend;
1400416022c9SBarry Smith   a->cstart     = oldmat->cstart;
1401416022c9SBarry Smith   a->cend       = oldmat->cend;
1402416022c9SBarry Smith   a->numtids    = oldmat->numtids;
1403416022c9SBarry Smith   a->mytid      = oldmat->mytid;
140464119d58SLois Curfman McInnes   a->insertmode = NOT_SET_VALUES;
1405d6dfbf8fSBarry Smith 
1406416022c9SBarry Smith   a->rowners    = (int *) PETSCMALLOC((a->numtids+1)*sizeof(int));CHKPTRQ(a->rowners);
1407416022c9SBarry Smith   PLogObjectMemory(mat,(a->numtids+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ));
1408416022c9SBarry Smith   PetscMemcpy(a->rowners,oldmat->rowners,(a->numtids+1)*sizeof(int));
1409416022c9SBarry Smith   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
14102ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
1411416022c9SBarry Smith     a->colmap      = (int *) PETSCMALLOC( (a->N)*sizeof(int) );CHKPTRQ(a->colmap);
1412416022c9SBarry Smith     PLogObjectMemory(mat,(a->N)*sizeof(int));
1413416022c9SBarry Smith     PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int));
1414416022c9SBarry Smith   } else a->colmap = 0;
1415ec8511deSBarry Smith   if (oldmat->garray && (len = ((Mat_SeqAIJ *) (oldmat->B->data))->n)) {
1416416022c9SBarry Smith     a->garray      = (int *) PETSCMALLOC(len*sizeof(int) ); CHKPTRQ(a->garray);
1417464493b3SBarry Smith     PLogObjectMemory(mat,len*sizeof(int));
1418416022c9SBarry Smith     PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int));
1419416022c9SBarry Smith   } else a->garray = 0;
1420d6dfbf8fSBarry Smith 
1421416022c9SBarry Smith   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
1422416022c9SBarry Smith   PLogObjectParent(mat,a->lvec);
1423416022c9SBarry Smith   ierr =  VecScatterCtxCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
1424416022c9SBarry Smith   PLogObjectParent(mat,a->Mvctx);
1425416022c9SBarry Smith   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
1426416022c9SBarry Smith   PLogObjectParent(mat,a->A);
1427416022c9SBarry Smith   ierr =  MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr);
1428416022c9SBarry Smith   PLogObjectParent(mat,a->B);
14298a729477SBarry Smith   *newmat = mat;
14308a729477SBarry Smith   return 0;
14318a729477SBarry Smith }
1432416022c9SBarry Smith 
1433416022c9SBarry Smith #include "sysio.h"
1434416022c9SBarry Smith 
143502834360SBarry Smith int MatLoad_MPIAIJorMPIRow(Viewer bview,MatType type,Mat *newmat)
1436416022c9SBarry Smith {
1437d65a2f8fSBarry Smith   Mat          A;
1438416022c9SBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
1439d65a2f8fSBarry Smith   Scalar       *vals,*svals;
1440416022c9SBarry Smith   PetscObject  vobj = (PetscObject) bview;
1441416022c9SBarry Smith   MPI_Comm     comm = vobj->comm;
1442416022c9SBarry Smith   MPI_Status   status;
1443416022c9SBarry Smith   int          header[4],mytid,numtid,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
1444d65a2f8fSBarry Smith   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
1445d65a2f8fSBarry Smith   int          tag = ((PetscObject)bview)->tag;
1446416022c9SBarry Smith 
1447416022c9SBarry Smith   MPI_Comm_size(comm,&numtid); MPI_Comm_rank(comm,&mytid);
1448416022c9SBarry Smith   if (!mytid) {
1449416022c9SBarry Smith     ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1450416022c9SBarry Smith     ierr = SYRead(fd,(char *)header,4,SYINT); CHKERRQ(ierr);
145102834360SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIAIJorMPIRow:not matrix object");
1452416022c9SBarry Smith   }
1453416022c9SBarry Smith 
1454416022c9SBarry Smith   MPI_Bcast(header+1,3,MPI_INT,0,comm);
1455416022c9SBarry Smith   M = header[1]; N = header[2];
1456416022c9SBarry Smith   /* determine ownership of all rows */
1457416022c9SBarry Smith   m = M/numtid + ((M % numtid) > mytid);
1458416022c9SBarry Smith   rowners = (int *) PETSCMALLOC((numtid+2)*sizeof(int)); CHKPTRQ(rowners);
1459416022c9SBarry Smith   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1460416022c9SBarry Smith   rowners[0] = 0;
1461416022c9SBarry Smith   for ( i=2; i<=numtid; i++ ) {
1462416022c9SBarry Smith     rowners[i] += rowners[i-1];
1463416022c9SBarry Smith   }
1464416022c9SBarry Smith   rstart = rowners[mytid];
1465416022c9SBarry Smith   rend   = rowners[mytid+1];
1466416022c9SBarry Smith 
1467416022c9SBarry Smith   /* distribute row lengths to all processors */
1468416022c9SBarry Smith   ourlens = (int*) PETSCMALLOC( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
1469416022c9SBarry Smith   offlens = ourlens + (rend-rstart);
1470416022c9SBarry Smith   if (!mytid) {
1471416022c9SBarry Smith     rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths);
1472416022c9SBarry Smith     ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
1473416022c9SBarry Smith     sndcounts = (int*) PETSCMALLOC( numtid*sizeof(int) ); CHKPTRQ(sndcounts);
1474416022c9SBarry Smith     for ( i=0; i<numtid; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1475d65a2f8fSBarry Smith     MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
14766ef5f50fSLois Curfman McInnes     PETSCFREE(sndcounts);
1477416022c9SBarry Smith   }
1478416022c9SBarry Smith   else {
1479416022c9SBarry Smith     MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);
1480416022c9SBarry Smith   }
1481416022c9SBarry Smith 
1482416022c9SBarry Smith   if (!mytid) {
1483416022c9SBarry Smith     /* calculate the number of nonzeros on each processor */
1484416022c9SBarry Smith     procsnz = (int*) PETSCMALLOC( numtid*sizeof(int) ); CHKPTRQ(procsnz);
1485416022c9SBarry Smith     PetscZero(procsnz,numtid*sizeof(int));
1486416022c9SBarry Smith     for ( i=0; i<numtid; i++ ) {
1487416022c9SBarry Smith       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
1488416022c9SBarry Smith         procsnz[i] += rowlengths[j];
1489416022c9SBarry Smith       }
1490416022c9SBarry Smith     }
14916ef5f50fSLois Curfman McInnes     PETSCFREE(rowlengths);
1492416022c9SBarry Smith 
1493416022c9SBarry Smith     /* determine max buffer needed and allocate it */
1494416022c9SBarry Smith     maxnz = 0;
1495416022c9SBarry Smith     for ( i=0; i<numtid; i++ ) {
1496416022c9SBarry Smith       maxnz = PETSCMAX(maxnz,procsnz[i]);
1497416022c9SBarry Smith     }
1498416022c9SBarry Smith     cols = (int *) PETSCMALLOC( maxnz*sizeof(int) ); CHKPTRQ(cols);
1499416022c9SBarry Smith 
1500416022c9SBarry Smith     /* read in my part of the matrix column indices  */
1501416022c9SBarry Smith     nz = procsnz[0];
1502d65a2f8fSBarry Smith     mycols = (int *) PETSCMALLOC( nz*sizeof(int) ); CHKPTRQ(mycols);
1503d65a2f8fSBarry Smith     ierr = SYRead(fd,mycols,nz,SYINT); CHKERRQ(ierr);
1504d65a2f8fSBarry Smith 
1505d65a2f8fSBarry Smith     /* read in every one elses and ship off */
1506d65a2f8fSBarry Smith     for ( i=1; i<numtid; i++ ) {
1507d65a2f8fSBarry Smith       nz = procsnz[i];
1508416022c9SBarry Smith       ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr);
1509d65a2f8fSBarry Smith       MPI_Send(cols,nz,MPI_INT,i,tag,comm);
1510d65a2f8fSBarry Smith     }
1511d65a2f8fSBarry Smith     PETSCFREE(cols);
1512416022c9SBarry Smith   }
1513416022c9SBarry Smith   else {
1514416022c9SBarry Smith     /* determine buffer space needed for message */
1515416022c9SBarry Smith     nz = 0;
1516416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
1517416022c9SBarry Smith       nz += ourlens[i];
1518416022c9SBarry Smith     }
1519d65a2f8fSBarry Smith     mycols = (int*) PETSCMALLOC( nz*sizeof(int) ); CHKPTRQ(mycols);
1520416022c9SBarry Smith 
1521416022c9SBarry Smith     /* receive message of column indices*/
1522d65a2f8fSBarry Smith     MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);
1523416022c9SBarry Smith     MPI_Get_count(&status,MPI_INT,&maxnz);
152402834360SBarry Smith     if (maxnz != nz) SETERRQ(1,"MatLoad_MPIAIJorMPIRow:something is wrong with file");
1525416022c9SBarry Smith   }
1526416022c9SBarry Smith 
1527416022c9SBarry Smith   /* loop over local rows, determining number of off diagonal entries */
1528416022c9SBarry Smith   PetscZero(offlens,m*sizeof(int));
1529416022c9SBarry Smith   jj = 0;
1530416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1531416022c9SBarry Smith     for ( j=0; j<ourlens[i]; j++ ) {
1532d65a2f8fSBarry Smith       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1533416022c9SBarry Smith       jj++;
1534416022c9SBarry Smith     }
1535416022c9SBarry Smith   }
1536d65a2f8fSBarry Smith 
1537d65a2f8fSBarry Smith 
1538d65a2f8fSBarry Smith   /* create our matrix */
1539416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1540416022c9SBarry Smith     ourlens[i] -= offlens[i];
1541416022c9SBarry Smith   }
1542d65a2f8fSBarry Smith   if (type == MATMPIAIJ) {
1543d65a2f8fSBarry Smith     ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr);
1544d65a2f8fSBarry Smith   }
1545d65a2f8fSBarry Smith   else if (type == MATMPIROW) {
1546d65a2f8fSBarry Smith     ierr = MatCreateMPIRow(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr);
1547d65a2f8fSBarry Smith   }
1548d65a2f8fSBarry Smith   A = *newmat;
1549d65a2f8fSBarry Smith   MatSetOption(A,COLUMNS_SORTED);
1550d65a2f8fSBarry Smith   for ( i=0; i<m; i++ ) {
1551d65a2f8fSBarry Smith     ourlens[i] += offlens[i];
1552d65a2f8fSBarry Smith   }
1553416022c9SBarry Smith 
1554416022c9SBarry Smith   if (!mytid) {
1555416022c9SBarry Smith     vals = (Scalar *) PETSCMALLOC( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
1556416022c9SBarry Smith 
1557416022c9SBarry Smith     /* read in my part of the matrix numerical values  */
1558416022c9SBarry Smith     nz = procsnz[0];
1559416022c9SBarry Smith     ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr);
1560d65a2f8fSBarry Smith 
1561d65a2f8fSBarry Smith     /* insert into matrix */
1562d65a2f8fSBarry Smith     jj      = rstart;
1563d65a2f8fSBarry Smith     smycols = mycols;
1564d65a2f8fSBarry Smith     svals   = vals;
1565d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1566d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1567d65a2f8fSBarry Smith       smycols += ourlens[i];
1568d65a2f8fSBarry Smith       svals   += ourlens[i];
1569d65a2f8fSBarry Smith       jj++;
1570416022c9SBarry Smith     }
1571416022c9SBarry Smith 
1572d65a2f8fSBarry Smith     /* read in other processors and ship out */
1573416022c9SBarry Smith     for ( i=1; i<numtid; i++ ) {
1574416022c9SBarry Smith       nz = procsnz[i];
1575416022c9SBarry Smith       ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr);
1576d65a2f8fSBarry Smith       MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
1577416022c9SBarry Smith     }
15786ef5f50fSLois Curfman McInnes     PETSCFREE(procsnz);
1579416022c9SBarry Smith   }
1580d65a2f8fSBarry Smith   else {
1581d65a2f8fSBarry Smith     /* receive numeric values */
1582d65a2f8fSBarry Smith     vals = (Scalar*) PETSCMALLOC( nz*sizeof(Scalar) ); CHKPTRQ(vals);
1583416022c9SBarry Smith 
1584d65a2f8fSBarry Smith     /* receive message of values*/
1585d65a2f8fSBarry Smith     MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
1586d65a2f8fSBarry Smith     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
158702834360SBarry Smith     if (maxnz != nz) SETERRQ(1,"MatLoad_MPIAIJorMPIRow:something is wrong with file");
1588d65a2f8fSBarry Smith 
1589d65a2f8fSBarry Smith     /* insert into matrix */
1590d65a2f8fSBarry Smith     jj      = rstart;
1591d65a2f8fSBarry Smith     smycols = mycols;
1592d65a2f8fSBarry Smith     svals   = vals;
1593d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1594d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1595d65a2f8fSBarry Smith       smycols += ourlens[i];
1596d65a2f8fSBarry Smith       svals   += ourlens[i];
1597d65a2f8fSBarry Smith       jj++;
1598d65a2f8fSBarry Smith     }
1599d65a2f8fSBarry Smith   }
16006ef5f50fSLois Curfman McInnes   PETSCFREE(ourlens); PETSCFREE(vals); PETSCFREE(mycols); PETSCFREE(rowners);
1601d65a2f8fSBarry Smith 
1602d65a2f8fSBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
1603d65a2f8fSBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
1604416022c9SBarry Smith   return 0;
1605416022c9SBarry Smith }
1606