xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision 9d00d63d4b4014b06eeb660baae4e576889960ed)
1cb512458SBarry Smith #ifndef lint
2*9d00d63dSBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.37 1995/05/03 01:03:31 curfman Exp bsmith $";
3cb512458SBarry Smith #endif
48a729477SBarry Smith 
51eb62cbbSBarry Smith #include "mpiaij.h"
68a729477SBarry Smith #include "vec/vecimpl.h"
7d6dfbf8fSBarry Smith #include "inline/spops.h"
88a729477SBarry Smith 
91eb62cbbSBarry Smith #define CHUNCKSIZE   100
101eb62cbbSBarry Smith /*
111eb62cbbSBarry Smith    This is a simple minded stash. Do a linear search to determine if
121eb62cbbSBarry Smith  in stash, if not add to end.
131eb62cbbSBarry Smith */
141eb62cbbSBarry Smith static int StashValues(Stash *stash,int row,int n, int *idxn,
151eb62cbbSBarry Smith                        Scalar *values,InsertMode addv)
168a729477SBarry Smith {
171eb62cbbSBarry Smith   int    i,j,N = stash->n,found,*n_idx, *n_idy;
181eb62cbbSBarry Smith   Scalar val,*n_array;
198a729477SBarry Smith 
201eb62cbbSBarry Smith   for ( i=0; i<n; i++ ) {
211eb62cbbSBarry Smith     found = 0;
221eb62cbbSBarry Smith     val = *values++;
238a729477SBarry Smith     for ( j=0; j<N; j++ ) {
241eb62cbbSBarry Smith       if ( stash->idx[j] == row && stash->idy[j] == idxn[i]) {
251eb62cbbSBarry Smith         /* found a match */
26*9d00d63dSBarry Smith         if (addv == ADDVALUES) stash->array[j] += val;
271eb62cbbSBarry Smith         else stash->array[j] = val;
281eb62cbbSBarry Smith         found = 1;
298a729477SBarry Smith         break;
308a729477SBarry Smith       }
318a729477SBarry Smith     }
321eb62cbbSBarry Smith     if (!found) { /* not found so add to end */
331eb62cbbSBarry Smith       if ( stash->n == stash->nmax ) {
341eb62cbbSBarry Smith         /* allocate a larger stash */
351eb62cbbSBarry Smith         n_array = (Scalar *) MALLOC( (stash->nmax + CHUNCKSIZE)*(
361eb62cbbSBarry Smith                                      2*sizeof(int) + sizeof(Scalar)));
371eb62cbbSBarry Smith         CHKPTR(n_array);
381eb62cbbSBarry Smith         n_idx = (int *) (n_array + stash->nmax + CHUNCKSIZE);
391eb62cbbSBarry Smith         n_idy = (int *) (n_idx + stash->nmax + CHUNCKSIZE);
401eb62cbbSBarry Smith         MEMCPY(n_array,stash->array,stash->nmax*sizeof(Scalar));
411eb62cbbSBarry Smith         MEMCPY(n_idx,stash->idx,stash->nmax*sizeof(int));
421eb62cbbSBarry Smith         MEMCPY(n_idy,stash->idy,stash->nmax*sizeof(int));
431eb62cbbSBarry Smith         if (stash->array) FREE(stash->array);
441eb62cbbSBarry Smith         stash->array = n_array; stash->idx = n_idx; stash->idy = n_idy;
451eb62cbbSBarry Smith         stash->nmax += CHUNCKSIZE;
461eb62cbbSBarry Smith       }
471eb62cbbSBarry Smith       stash->array[stash->n]   = val;
481eb62cbbSBarry Smith       stash->idx[stash->n]     = row;
491eb62cbbSBarry Smith       stash->idy[stash->n++]   = idxn[i];
501eb62cbbSBarry Smith     }
518a729477SBarry Smith   }
528a729477SBarry Smith   return 0;
538a729477SBarry Smith }
548a729477SBarry Smith 
559e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column
569e25ed09SBarry Smith number to the local number in the off-diagonal part of the local
579e25ed09SBarry Smith storage of the matrix.  This is done in a non scable way since the
589e25ed09SBarry Smith length of colmap equals the global matrix length.
599e25ed09SBarry Smith */
609e25ed09SBarry Smith static int CreateColmap(Mat mat)
619e25ed09SBarry Smith {
6244a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
6344a69424SLois Curfman McInnes   Mat_AIJ    *B = (Mat_AIJ*) aij->B->data;
649e25ed09SBarry Smith   int        n = B->n,i;
659e25ed09SBarry Smith   aij->colmap = (int *) MALLOC( aij->N*sizeof(int) ); CHKPTR(aij->colmap);
669e25ed09SBarry Smith   MEMSET(aij->colmap,0,aij->N*sizeof(int));
67abc0e9e4SLois Curfman McInnes   for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1;
689e25ed09SBarry Smith   return 0;
699e25ed09SBarry Smith }
709e25ed09SBarry Smith 
712ee70a88SLois Curfman McInnes static int MatSetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,
721eb62cbbSBarry Smith                             int *idxn,Scalar *v,InsertMode addv)
738a729477SBarry Smith {
7444a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
751eb62cbbSBarry Smith   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
761eb62cbbSBarry Smith   int        cstart = aij->cstart, cend = aij->cend,row,col;
778a729477SBarry Smith 
7807a0e7adSLois Curfman McInnes   if (aij->insertmode != NOTSETVALUES && aij->insertmode != addv) {
791eb62cbbSBarry Smith     SETERR(1,"You cannot mix inserts and adds");
808a729477SBarry Smith   }
811eb62cbbSBarry Smith   aij->insertmode = addv;
828a729477SBarry Smith   for ( i=0; i<m; i++ ) {
83da3a660dSBarry Smith     if (idxm[i] < 0) SETERR(1,"Negative row index");
84da3a660dSBarry Smith     if (idxm[i] >= aij->M) SETERR(1,"Row index too large");
851eb62cbbSBarry Smith     if (idxm[i] >= rstart && idxm[i] < rend) {
861eb62cbbSBarry Smith       row = idxm[i] - rstart;
871eb62cbbSBarry Smith       for ( j=0; j<n; j++ ) {
88da3a660dSBarry Smith         if (idxn[j] < 0) SETERR(1,"Negative column index");
89da3a660dSBarry Smith         if (idxn[j] >= aij->N) SETERR(1,"Column index too large");
901eb62cbbSBarry Smith         if (idxn[j] >= cstart && idxn[j] < cend){
911eb62cbbSBarry Smith           col = idxn[j] - cstart;
921eb62cbbSBarry Smith           ierr = MatSetValues(aij->A,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr);
931eb62cbbSBarry Smith         }
941eb62cbbSBarry Smith         else {
95d6dfbf8fSBarry Smith           if (aij->assembled) {
969e25ed09SBarry Smith             if (!aij->colmap) {ierr = CreateColmap(mat); CHKERR(ierr);}
979e25ed09SBarry Smith             col = aij->colmap[idxn[j]] - 1;
989e25ed09SBarry Smith             if (col < 0) {
999e25ed09SBarry Smith               SETERR(1,"Cannot insert new off diagonal block nonzero in\
1009e25ed09SBarry Smith                      already\
101d6dfbf8fSBarry Smith                      assembled matrix. Contact petsc-maint@mcs.anl.gov\
102d6dfbf8fSBarry Smith                      if your need this feature");
103d6dfbf8fSBarry Smith             }
1049e25ed09SBarry Smith           }
1059e25ed09SBarry Smith           else col = idxn[j];
1061eb62cbbSBarry Smith           ierr = MatSetValues(aij->B,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr);
1071eb62cbbSBarry Smith         }
1081eb62cbbSBarry Smith       }
1091eb62cbbSBarry Smith     }
1101eb62cbbSBarry Smith     else {
1111eb62cbbSBarry Smith       ierr = StashValues(&aij->stash,idxm[i],n,idxn,v+i*n,addv);CHKERR(ierr);
1121eb62cbbSBarry Smith     }
1138a729477SBarry Smith   }
1148a729477SBarry Smith   return 0;
1158a729477SBarry Smith }
1168a729477SBarry Smith 
1178a729477SBarry Smith /*
1181eb62cbbSBarry Smith     the assembly code is alot like the code for vectors, we should
1191eb62cbbSBarry Smith     sometime derive a single assembly code that can be used for
1201eb62cbbSBarry Smith     either case.
1218a729477SBarry Smith */
1228a729477SBarry Smith 
123fc76ce87SLois Curfman McInnes static int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode)
1248a729477SBarry Smith {
12544a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
126d6dfbf8fSBarry Smith   MPI_Comm    comm = mat->comm;
1276abc6512SBarry Smith   int         numtids = aij->numtids, *owners = aij->rowners;
1281eb62cbbSBarry Smith   int         mytid = aij->mytid;
1291eb62cbbSBarry Smith   MPI_Request *send_waits,*recv_waits;
1306abc6512SBarry Smith   int         *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
1311eb62cbbSBarry Smith   int         tag = 50, *owner,*starts,count;
1321eb62cbbSBarry Smith   InsertMode  addv;
1331eb62cbbSBarry Smith   Scalar      *rvalues,*svalues;
1341eb62cbbSBarry Smith 
1351eb62cbbSBarry Smith   /* make sure all processors are either in INSERTMODE or ADDMODE */
13628988994SBarry Smith   MPI_Allreduce((void *) &aij->insertmode,(void *) &addv,1,MPI_INT,
1371eb62cbbSBarry Smith                 MPI_BOR,comm);
138*9d00d63dSBarry Smith   if (addv == (ADDVALUES|INSERTVALUES)) {
1391eb62cbbSBarry Smith     SETERR(1,"Some processors have inserted while others have added");
1401eb62cbbSBarry Smith   }
1411eb62cbbSBarry Smith   aij->insertmode = addv; /* in case this processor had no cache */
1421eb62cbbSBarry Smith 
1431eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
1441eb62cbbSBarry Smith   nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs);
1451eb62cbbSBarry Smith   MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids;
1461eb62cbbSBarry Smith   owner = (int *) MALLOC( (aij->stash.n+1)*sizeof(int) ); CHKPTR(owner);
1471eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
1481eb62cbbSBarry Smith     idx = aij->stash.idx[i];
1491eb62cbbSBarry Smith     for ( j=0; j<numtids; j++ ) {
1501eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
1511eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
1528a729477SBarry Smith       }
1538a729477SBarry Smith     }
1548a729477SBarry Smith   }
1551eb62cbbSBarry Smith   nsends = 0;  for ( i=0; i<numtids; i++ ) { nsends += procs[i];}
1561eb62cbbSBarry Smith 
1571eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
1581eb62cbbSBarry Smith   work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work);
1591eb62cbbSBarry Smith   MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm);
1601eb62cbbSBarry Smith   nreceives = work[mytid];
1611eb62cbbSBarry Smith   MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm);
1621eb62cbbSBarry Smith   nmax = work[mytid];
1631eb62cbbSBarry Smith   FREE(work);
1641eb62cbbSBarry Smith 
1651eb62cbbSBarry Smith   /* post receives:
1661eb62cbbSBarry Smith        1) each message will consist of ordered pairs
1671eb62cbbSBarry Smith      (global index,value) we store the global index as a double
168d6dfbf8fSBarry Smith      to simplify the message passing.
1691eb62cbbSBarry Smith        2) since we don't know how long each individual message is we
1701eb62cbbSBarry Smith      allocate the largest needed buffer for each receive. Potentially
1711eb62cbbSBarry Smith      this is a lot of wasted space.
1721eb62cbbSBarry Smith 
1731eb62cbbSBarry Smith 
1741eb62cbbSBarry Smith        This could be done better.
1751eb62cbbSBarry Smith   */
17628988994SBarry Smith   rvalues = (Scalar *) MALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));
1771eb62cbbSBarry Smith   CHKPTR(rvalues);
1781eb62cbbSBarry Smith   recv_waits = (MPI_Request *) MALLOC((nreceives+1)*sizeof(MPI_Request));
1791eb62cbbSBarry Smith   CHKPTR(recv_waits);
1801eb62cbbSBarry Smith   for ( i=0; i<nreceives; i++ ) {
1811eb62cbbSBarry Smith     MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPI_SCALAR,MPI_ANY_SOURCE,tag,
1821eb62cbbSBarry Smith               comm,recv_waits+i);
1831eb62cbbSBarry Smith   }
1841eb62cbbSBarry Smith 
1851eb62cbbSBarry Smith   /* do sends:
1861eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
1871eb62cbbSBarry Smith          the ith processor
1881eb62cbbSBarry Smith   */
1891eb62cbbSBarry Smith   svalues = (Scalar *) MALLOC( 3*(aij->stash.n+1)*sizeof(Scalar) );
1901eb62cbbSBarry Smith   CHKPTR(svalues);
1911eb62cbbSBarry Smith   send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request));
1921eb62cbbSBarry Smith   CHKPTR(send_waits);
1931eb62cbbSBarry Smith   starts = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(starts);
1941eb62cbbSBarry Smith   starts[0] = 0;
1951eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1961eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
1971eb62cbbSBarry Smith     svalues[3*starts[owner[i]]]       = (Scalar)  aij->stash.idx[i];
1981eb62cbbSBarry Smith     svalues[3*starts[owner[i]]+1]     = (Scalar)  aij->stash.idy[i];
1991eb62cbbSBarry Smith     svalues[3*(starts[owner[i]]++)+2] =  aij->stash.array[i];
2001eb62cbbSBarry Smith   }
2011eb62cbbSBarry Smith   FREE(owner);
2021eb62cbbSBarry Smith   starts[0] = 0;
2031eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
2041eb62cbbSBarry Smith   count = 0;
2051eb62cbbSBarry Smith   for ( i=0; i<numtids; i++ ) {
2061eb62cbbSBarry Smith     if (procs[i]) {
2071eb62cbbSBarry Smith       MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPI_SCALAR,i,tag,
2081eb62cbbSBarry Smith                 comm,send_waits+count++);
2091eb62cbbSBarry Smith     }
2101eb62cbbSBarry Smith   }
2111eb62cbbSBarry Smith   FREE(starts); FREE(nprocs);
2121eb62cbbSBarry Smith 
2131eb62cbbSBarry Smith   /* Free cache space */
2141eb62cbbSBarry Smith   aij->stash.nmax = aij->stash.n = 0;
2151eb62cbbSBarry Smith   if (aij->stash.array){ FREE(aij->stash.array); aij->stash.array = 0;}
2161eb62cbbSBarry Smith 
2171eb62cbbSBarry Smith   aij->svalues    = svalues;       aij->rvalues = rvalues;
2181eb62cbbSBarry Smith   aij->nsends     = nsends;         aij->nrecvs = nreceives;
2191eb62cbbSBarry Smith   aij->send_waits = send_waits; aij->recv_waits = recv_waits;
2201eb62cbbSBarry Smith   aij->rmax       = nmax;
2211eb62cbbSBarry Smith 
2221eb62cbbSBarry Smith   return 0;
2231eb62cbbSBarry Smith }
22444a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat);
2251eb62cbbSBarry Smith 
226fc76ce87SLois Curfman McInnes static int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode)
2271eb62cbbSBarry Smith {
2281eb62cbbSBarry Smith   int        ierr;
22944a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
2301eb62cbbSBarry Smith 
2311eb62cbbSBarry Smith   MPI_Status  *send_status,recv_status;
2326abc6512SBarry Smith   int         imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n;
2331eb62cbbSBarry Smith   int         row,col;
2341eb62cbbSBarry Smith   Scalar      *values,val;
2351eb62cbbSBarry Smith   InsertMode  addv = aij->insertmode;
2361eb62cbbSBarry Smith 
2371eb62cbbSBarry Smith   /*  wait on receives */
2381eb62cbbSBarry Smith   while (count) {
239d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);
2401eb62cbbSBarry Smith     /* unpack receives into our local space */
241d6dfbf8fSBarry Smith     values = aij->rvalues + 3*imdex*aij->rmax;
2421eb62cbbSBarry Smith     MPI_Get_count(&recv_status,MPI_SCALAR,&n);
2431eb62cbbSBarry Smith     n = n/3;
2441eb62cbbSBarry Smith     for ( i=0; i<n; i++ ) {
2451eb62cbbSBarry Smith       row = (int) PETSCREAL(values[3*i]) - aij->rstart;
2461eb62cbbSBarry Smith       col = (int) PETSCREAL(values[3*i+1]);
2471eb62cbbSBarry Smith       val = values[3*i+2];
2481eb62cbbSBarry Smith       if (col >= aij->cstart && col < aij->cend) {
2491eb62cbbSBarry Smith           col -= aij->cstart;
2501eb62cbbSBarry Smith         MatSetValues(aij->A,1,&row,1,&col,&val,addv);
2511eb62cbbSBarry Smith       }
2521eb62cbbSBarry Smith       else {
253d6dfbf8fSBarry Smith         if (aij->assembled) {
2549e25ed09SBarry Smith           if (!aij->colmap) {ierr = CreateColmap(mat); CHKERR(ierr);}
2559e25ed09SBarry Smith           col = aij->colmap[col] - 1;
2569e25ed09SBarry Smith           if (col < 0) {
2579e25ed09SBarry Smith             SETERR(1,"Cannot insert new off diagonal block nonzero in\
2589e25ed09SBarry Smith                      already\
259d6dfbf8fSBarry Smith                      assembled matrix. Contact petsc-maint@mcs.anl.gov\
260d6dfbf8fSBarry Smith                      if your need this feature");
261d6dfbf8fSBarry Smith           }
2629e25ed09SBarry Smith         }
2631eb62cbbSBarry Smith         MatSetValues(aij->B,1,&row,1,&col,&val,addv);
2641eb62cbbSBarry Smith       }
2651eb62cbbSBarry Smith     }
2661eb62cbbSBarry Smith     count--;
2671eb62cbbSBarry Smith   }
2681eb62cbbSBarry Smith   FREE(aij->recv_waits); FREE(aij->rvalues);
2691eb62cbbSBarry Smith 
2701eb62cbbSBarry Smith   /* wait on sends */
2711eb62cbbSBarry Smith   if (aij->nsends) {
2721eb62cbbSBarry Smith     send_status = (MPI_Status *) MALLOC( aij->nsends*sizeof(MPI_Status) );
2731eb62cbbSBarry Smith     CHKPTR(send_status);
2741eb62cbbSBarry Smith     MPI_Waitall(aij->nsends,aij->send_waits,send_status);
2751eb62cbbSBarry Smith     FREE(send_status);
2761eb62cbbSBarry Smith   }
2771eb62cbbSBarry Smith   FREE(aij->send_waits); FREE(aij->svalues);
2781eb62cbbSBarry Smith 
27907a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
280ee50ffe9SBarry Smith   ierr = MatAssemblyBegin(aij->A,mode); CHKERR(ierr);
281ee50ffe9SBarry Smith   ierr = MatAssemblyEnd(aij->A,mode); CHKERR(ierr);
2821eb62cbbSBarry Smith 
2835e42470aSBarry Smith   if (!aij->assembled && mode == FINAL_ASSEMBLY) {
284d3c9fed9SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERR(ierr);
285d6dfbf8fSBarry Smith     aij->assembled = 1;
2865e42470aSBarry Smith   }
287ee50ffe9SBarry Smith   ierr = MatAssemblyBegin(aij->B,mode); CHKERR(ierr);
288ee50ffe9SBarry Smith   ierr = MatAssemblyEnd(aij->B,mode); CHKERR(ierr);
2895e42470aSBarry Smith 
2908a729477SBarry Smith   return 0;
2918a729477SBarry Smith }
2928a729477SBarry Smith 
2932ee70a88SLois Curfman McInnes static int MatZeroEntries_MPIAIJ(Mat A)
2941eb62cbbSBarry Smith {
29544a69424SLois Curfman McInnes   Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data;
2961eb62cbbSBarry Smith 
2971eb62cbbSBarry Smith   MatZeroEntries(l->A); MatZeroEntries(l->B);
2981eb62cbbSBarry Smith   return 0;
2991eb62cbbSBarry Smith }
3001eb62cbbSBarry Smith 
3011eb62cbbSBarry Smith /* again this uses the same basic stratagy as in the assembly and
3021eb62cbbSBarry Smith    scatter create routines, we should try to do it systemamatically
3031eb62cbbSBarry Smith    if we can figure out the proper level of generality. */
3041eb62cbbSBarry Smith 
3051eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the
3061eb62cbbSBarry Smith    matrix is square and the column and row owerships are identical.
3071eb62cbbSBarry Smith    This is a BUG. The only way to fix it seems to be to access
3081eb62cbbSBarry Smith    aij->A and aij->B directly and not through the MatZeroRows()
3091eb62cbbSBarry Smith    routine.
3101eb62cbbSBarry Smith */
31144a69424SLois Curfman McInnes static int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag)
3121eb62cbbSBarry Smith {
31344a69424SLois Curfman McInnes   Mat_MPIAIJ     *l = (Mat_MPIAIJ *) A->data;
3141eb62cbbSBarry Smith   int            i,ierr,N, *rows,*owners = l->rowners,numtids = l->numtids;
3156abc6512SBarry Smith   int            *procs,*nprocs,j,found,idx,nsends,*work;
3161eb62cbbSBarry Smith   int            nmax,*svalues,*starts,*owner,nrecvs,mytid = l->mytid;
3176abc6512SBarry Smith   int            *rvalues,tag = 67,count,base,slen,n,*source;
318d6dfbf8fSBarry Smith   int            *lens,imdex,*lrows,*values;
319d6dfbf8fSBarry Smith   MPI_Comm       comm = A->comm;
3201eb62cbbSBarry Smith   MPI_Request    *send_waits,*recv_waits;
3211eb62cbbSBarry Smith   MPI_Status     recv_status,*send_status;
3221eb62cbbSBarry Smith   IS             istmp;
3231eb62cbbSBarry Smith 
32444a69424SLois Curfman McInnes   if (!l->assembled) SETERR(1,"MatZeroRows_MPIAIJ: must assemble matrix first");
3251eb62cbbSBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERR(ierr);
3261eb62cbbSBarry Smith   ierr = ISGetIndices(is,&rows); CHKERR(ierr);
3271eb62cbbSBarry Smith 
3281eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
3291eb62cbbSBarry Smith   nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs);
3301eb62cbbSBarry Smith   MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids;
3311eb62cbbSBarry Smith   owner = (int *) MALLOC((N+1)*sizeof(int)); CHKPTR(owner); /* see note*/
3321eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
3331eb62cbbSBarry Smith     idx = rows[i];
3341eb62cbbSBarry Smith     found = 0;
3351eb62cbbSBarry Smith     for ( j=0; j<numtids; j++ ) {
3361eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
3371eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
3381eb62cbbSBarry Smith       }
3391eb62cbbSBarry Smith     }
340d6dfbf8fSBarry Smith     if (!found) SETERR(1,"Imdex out of range");
3411eb62cbbSBarry Smith   }
3421eb62cbbSBarry Smith   nsends = 0;  for ( i=0; i<numtids; i++ ) { nsends += procs[i];}
3431eb62cbbSBarry Smith 
3441eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
3451eb62cbbSBarry Smith   work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work);
3461eb62cbbSBarry Smith   MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm);
3471eb62cbbSBarry Smith   nrecvs = work[mytid];
3481eb62cbbSBarry Smith   MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm);
3491eb62cbbSBarry Smith   nmax = work[mytid];
3501eb62cbbSBarry Smith   FREE(work);
3511eb62cbbSBarry Smith 
3521eb62cbbSBarry Smith   /* post receives:   */
35328988994SBarry Smith   rvalues = (int *) MALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
3541eb62cbbSBarry Smith   CHKPTR(rvalues);
3551eb62cbbSBarry Smith   recv_waits = (MPI_Request *) MALLOC((nrecvs+1)*sizeof(MPI_Request));
3561eb62cbbSBarry Smith   CHKPTR(recv_waits);
3571eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
3581eb62cbbSBarry Smith     MPI_Irecv((void *)(rvalues+nmax*i),nmax,MPI_INT,MPI_ANY_SOURCE,tag,
3591eb62cbbSBarry Smith               comm,recv_waits+i);
3601eb62cbbSBarry Smith   }
3611eb62cbbSBarry Smith 
3621eb62cbbSBarry Smith   /* do sends:
3631eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
3641eb62cbbSBarry Smith          the ith processor
3651eb62cbbSBarry Smith   */
3661eb62cbbSBarry Smith   svalues = (int *) MALLOC( (N+1)*sizeof(int) ); CHKPTR(svalues);
3671eb62cbbSBarry Smith   send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request));
3681eb62cbbSBarry Smith   CHKPTR(send_waits);
3691eb62cbbSBarry Smith   starts = (int *) MALLOC( (numtids+1)*sizeof(int) ); CHKPTR(starts);
3701eb62cbbSBarry Smith   starts[0] = 0;
3711eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3721eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
3731eb62cbbSBarry Smith     svalues[starts[owner[i]]++] = rows[i];
3741eb62cbbSBarry Smith   }
3751eb62cbbSBarry Smith   ISRestoreIndices(is,&rows);
3761eb62cbbSBarry Smith 
3771eb62cbbSBarry Smith   starts[0] = 0;
3781eb62cbbSBarry Smith   for ( i=1; i<numtids+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3791eb62cbbSBarry Smith   count = 0;
3801eb62cbbSBarry Smith   for ( i=0; i<numtids; i++ ) {
3811eb62cbbSBarry Smith     if (procs[i]) {
3821eb62cbbSBarry Smith       MPI_Isend((void*)(svalues+starts[i]),nprocs[i],MPI_INT,i,tag,
3831eb62cbbSBarry Smith                 comm,send_waits+count++);
3841eb62cbbSBarry Smith     }
3851eb62cbbSBarry Smith   }
3861eb62cbbSBarry Smith   FREE(starts);
3871eb62cbbSBarry Smith 
3881eb62cbbSBarry Smith   base = owners[mytid];
3891eb62cbbSBarry Smith 
3901eb62cbbSBarry Smith   /*  wait on receives */
3911eb62cbbSBarry Smith   lens = (int *) MALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTR(lens);
3921eb62cbbSBarry Smith   source = lens + nrecvs;
3931eb62cbbSBarry Smith   count = nrecvs; slen = 0;
3941eb62cbbSBarry Smith   while (count) {
395d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
3961eb62cbbSBarry Smith     /* unpack receives into our local space */
3971eb62cbbSBarry Smith     MPI_Get_count(&recv_status,MPI_INT,&n);
398d6dfbf8fSBarry Smith     source[imdex]  = recv_status.MPI_SOURCE;
399d6dfbf8fSBarry Smith     lens[imdex]  = n;
4001eb62cbbSBarry Smith     slen += n;
4011eb62cbbSBarry Smith     count--;
4021eb62cbbSBarry Smith   }
4031eb62cbbSBarry Smith   FREE(recv_waits);
4041eb62cbbSBarry Smith 
4051eb62cbbSBarry Smith   /* move the data into the send scatter */
4061eb62cbbSBarry Smith   lrows = (int *) MALLOC( slen*sizeof(int) ); CHKPTR(lrows);
4071eb62cbbSBarry Smith   count = 0;
4081eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
4091eb62cbbSBarry Smith     values = rvalues + i*nmax;
4101eb62cbbSBarry Smith     for ( j=0; j<lens[i]; j++ ) {
4111eb62cbbSBarry Smith       lrows[count++] = values[j] - base;
4121eb62cbbSBarry Smith     }
4131eb62cbbSBarry Smith   }
4141eb62cbbSBarry Smith   FREE(rvalues); FREE(lens);
4151eb62cbbSBarry Smith   FREE(owner); FREE(nprocs);
4161eb62cbbSBarry Smith 
4171eb62cbbSBarry Smith   /* actually zap the local rows */
4186b5873e3SBarry Smith   ierr = ISCreateSequential(MPI_COMM_SELF,slen,lrows,&istmp);
4196b5873e3SBarry Smith   CHKERR(ierr);  FREE(lrows);
4201eb62cbbSBarry Smith   ierr = MatZeroRows(l->A,istmp,diag); CHKERR(ierr);
4211eb62cbbSBarry Smith   ierr = MatZeroRows(l->B,istmp,0); CHKERR(ierr);
4221eb62cbbSBarry Smith   ierr = ISDestroy(istmp); CHKERR(ierr);
4231eb62cbbSBarry Smith 
4241eb62cbbSBarry Smith   /* wait on sends */
4251eb62cbbSBarry Smith   if (nsends) {
4261eb62cbbSBarry Smith     send_status = (MPI_Status *) MALLOC( nsends*sizeof(MPI_Status) );
4271eb62cbbSBarry Smith     CHKPTR(send_status);
4281eb62cbbSBarry Smith     MPI_Waitall(nsends,send_waits,send_status);
4291eb62cbbSBarry Smith     FREE(send_status);
4301eb62cbbSBarry Smith   }
4311eb62cbbSBarry Smith   FREE(send_waits); FREE(svalues);
4321eb62cbbSBarry Smith 
4331eb62cbbSBarry Smith 
4341eb62cbbSBarry Smith   return 0;
4351eb62cbbSBarry Smith }
4361eb62cbbSBarry Smith 
43744a69424SLois Curfman McInnes static int MatMult_MPIAIJ(Mat aijin,Vec xx,Vec yy)
4381eb62cbbSBarry Smith {
43944a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
4401eb62cbbSBarry Smith   int        ierr;
44144a69424SLois Curfman McInnes   if (!aij->assembled) SETERR(1,"MatMult_MPIAIJ: must assemble matrix first");
442*9d00d63dSBarry Smith   ierr = VecScatterBegin(xx,0,aij->lvec,0,INSERTVALUES,SCATTERALL,aij->Mvctx);
4431eb62cbbSBarry Smith   CHKERR(ierr);
4441eb62cbbSBarry Smith   ierr = MatMult(aij->A,xx,yy); CHKERR(ierr);
445*9d00d63dSBarry Smith   ierr = VecScatterEnd(xx,0,aij->lvec,0,INSERTVALUES,SCATTERALL,aij->Mvctx);
446d6dfbf8fSBarry Smith   CHKERR(ierr);
4471eb62cbbSBarry Smith   ierr = MatMultAdd(aij->B,aij->lvec,yy,yy); CHKERR(ierr);
4481eb62cbbSBarry Smith   return 0;
4491eb62cbbSBarry Smith }
4501eb62cbbSBarry Smith 
45144a69424SLois Curfman McInnes static int MatMultAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz)
452da3a660dSBarry Smith {
45344a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
454da3a660dSBarry Smith   int        ierr;
45544a69424SLois Curfman McInnes   if (!aij->assembled) SETERR(1,"MatMult_MPIAIJ: must assemble matrix first");
456*9d00d63dSBarry Smith   ierr = VecScatterBegin(xx,0,aij->lvec,0,INSERTVALUES,SCATTERALL,aij->Mvctx);
457da3a660dSBarry Smith   CHKERR(ierr);
458da3a660dSBarry Smith   ierr = MatMultAdd(aij->A,xx,yy,zz); CHKERR(ierr);
459*9d00d63dSBarry Smith   ierr = VecScatterEnd(xx,0,aij->lvec,0,INSERTVALUES,SCATTERALL,aij->Mvctx);
460da3a660dSBarry Smith   CHKERR(ierr);
461da3a660dSBarry Smith   ierr = MatMultAdd(aij->B,aij->lvec,zz,zz); CHKERR(ierr);
462da3a660dSBarry Smith   return 0;
463da3a660dSBarry Smith }
464da3a660dSBarry Smith 
46544a69424SLois Curfman McInnes static int MatMultTrans_MPIAIJ(Mat aijin,Vec xx,Vec yy)
466da3a660dSBarry Smith {
46744a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
468da3a660dSBarry Smith   int        ierr;
469da3a660dSBarry Smith 
47044a69424SLois Curfman McInnes   if (!aij->assembled)
47144a69424SLois Curfman McInnes     SETERR(1,"MatMulTrans_MPIAIJ: must assemble matrix first");
472da3a660dSBarry Smith   /* do nondiagonal part */
473da3a660dSBarry Smith   ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr);
474da3a660dSBarry Smith   /* send it on its way */
475*9d00d63dSBarry Smith   ierr = VecScatterBegin(aij->lvec,0,yy,0,ADDVALUES,
476*9d00d63dSBarry Smith                          SCATTERALL|SCATTERREVERSE,aij->Mvctx); CHKERR(ierr);
477da3a660dSBarry Smith   /* do local part */
478da3a660dSBarry Smith   ierr = MatMultTrans(aij->A,xx,yy); CHKERR(ierr);
479da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
480da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
481da3a660dSBarry Smith   /* but is not perhaps always true. */
482*9d00d63dSBarry Smith   ierr = VecScatterEnd(aij->lvec,0,yy,0,ADDVALUES,SCATTERALL|SCATTERREVERSE,
483da3a660dSBarry Smith                          aij->Mvctx); CHKERR(ierr);
484da3a660dSBarry Smith   return 0;
485da3a660dSBarry Smith }
486da3a660dSBarry Smith 
48744a69424SLois Curfman McInnes static int MatMultTransAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz)
488da3a660dSBarry Smith {
48944a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
490da3a660dSBarry Smith   int        ierr;
491da3a660dSBarry Smith 
49244a69424SLois Curfman McInnes   if (!aij->assembled)
49344a69424SLois Curfman McInnes     SETERR(1,"MatMulTransAdd_MPIAIJ: must assemble matrix first");
494da3a660dSBarry Smith   /* do nondiagonal part */
495da3a660dSBarry Smith   ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr);
496da3a660dSBarry Smith   /* send it on its way */
497*9d00d63dSBarry Smith   ierr = VecScatterBegin(aij->lvec,0,zz,0,ADDVALUES,
498*9d00d63dSBarry Smith                          SCATTERALL|SCATTERREVERSE,aij->Mvctx); CHKERR(ierr);
499da3a660dSBarry Smith   /* do local part */
500da3a660dSBarry Smith   ierr = MatMultTransAdd(aij->A,xx,yy,zz); CHKERR(ierr);
501da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
502da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
503da3a660dSBarry Smith   /* but is not perhaps always true. */
504*9d00d63dSBarry Smith   ierr = VecScatterEnd(aij->lvec,0,zz,0,ADDVALUES,SCATTERALL|SCATTERREVERSE,
505da3a660dSBarry Smith                          aij->Mvctx); CHKERR(ierr);
506da3a660dSBarry Smith   return 0;
507da3a660dSBarry Smith }
508da3a660dSBarry Smith 
5091eb62cbbSBarry Smith /*
5101eb62cbbSBarry Smith   This only works correctly for square matrices where the subblock A->A is the
5111eb62cbbSBarry Smith    diagonal block
5121eb62cbbSBarry Smith */
513d1710a03SLois Curfman McInnes static int MatGetDiagonal_MPIAIJ(Mat Ain,Vec v)
5141eb62cbbSBarry Smith {
51544a69424SLois Curfman McInnes   Mat_MPIAIJ *A = (Mat_MPIAIJ *) Ain->data;
51644a69424SLois Curfman McInnes   if (!A->assembled) SETERR(1,"MatGetDiag_MPIAIJ: must assemble matrix first");
5171eb62cbbSBarry Smith   return MatGetDiagonal(A->A,v);
5181eb62cbbSBarry Smith }
5191eb62cbbSBarry Smith 
52044a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj)
5211eb62cbbSBarry Smith {
5221eb62cbbSBarry Smith   Mat        mat = (Mat) obj;
52344a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
5241eb62cbbSBarry Smith   int        ierr;
525a5a9c739SBarry Smith #if defined(PETSC_LOG)
526a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",aij->M,aij->N);
527a5a9c739SBarry Smith #endif
5281eb62cbbSBarry Smith   FREE(aij->rowners);
5291eb62cbbSBarry Smith   ierr = MatDestroy(aij->A); CHKERR(ierr);
5301eb62cbbSBarry Smith   ierr = MatDestroy(aij->B); CHKERR(ierr);
5319e25ed09SBarry Smith   if (aij->colmap) FREE(aij->colmap);
5329e25ed09SBarry Smith   if (aij->garray) FREE(aij->garray);
5331eb62cbbSBarry Smith   if (aij->lvec) VecDestroy(aij->lvec);
5341eb62cbbSBarry Smith   if (aij->Mvctx) VecScatterCtxDestroy(aij->Mvctx);
535a5a9c739SBarry Smith   FREE(aij);
536a5a9c739SBarry Smith   PLogObjectDestroy(mat);
537a5a9c739SBarry Smith   PETSCHEADERDESTROY(mat);
5381eb62cbbSBarry Smith   return 0;
5391eb62cbbSBarry Smith }
5407857610eSBarry Smith #include "draw.h"
541ee50ffe9SBarry Smith #include "pviewer.h"
542ee50ffe9SBarry Smith 
54344a69424SLois Curfman McInnes static int MatView_MPIAIJ(PetscObject obj,Viewer viewer)
5441eb62cbbSBarry Smith {
5451eb62cbbSBarry Smith   Mat        mat = (Mat) obj;
54644a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
5471eb62cbbSBarry Smith   int        ierr;
548d13ab20cSBarry Smith   PetscObject vobj = (PetscObject) viewer;
5491eb62cbbSBarry Smith 
550ab254492SBarry Smith   if (!aij->assembled) SETERR(1,"MatView_MPIAIJ: must assemble matrix first");
551154123eaSLois Curfman McInnes   if (!viewer) { /* so that viewers may be used from debuggers */
552154123eaSLois Curfman McInnes     viewer = STDOUT_VIEWER; vobj = (PetscObject) viewer;
553154123eaSLois Curfman McInnes   }
554ab254492SBarry Smith   if (vobj->cookie == DRAW_COOKIE && vobj->type == NULLWINDOW) return 0;
5557857610eSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == FILE_VIEWER) {
5564a0ce102SLois Curfman McInnes     FILE *fd = ViewerFileGetPointer_Private(viewer);
557d6dfbf8fSBarry Smith     MPE_Seq_begin(mat->comm,1);
558d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
5591eb62cbbSBarry Smith              aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
5601eb62cbbSBarry Smith              aij->cend);
561da69df5fSBarry Smith     ierr = MatView(aij->A,viewer); CHKERR(ierr);
562da69df5fSBarry Smith     ierr = MatView(aij->B,viewer); CHKERR(ierr);
563d13ab20cSBarry Smith     fflush(fd);
564d6dfbf8fSBarry Smith     MPE_Seq_end(mat->comm,1);
565d13ab20cSBarry Smith   }
5667857610eSBarry Smith   else if ((vobj->cookie == VIEWER_COOKIE && vobj->type == FILES_VIEWER) ||
5677857610eSBarry Smith             vobj->cookie == DRAW_COOKIE) {
56895373324SBarry Smith     int numtids = aij->numtids, mytid = aij->mytid;
56995373324SBarry Smith     if (numtids == 1) {
57095373324SBarry Smith       ierr = MatView(aij->A,viewer); CHKERR(ierr);
57195373324SBarry Smith     }
57295373324SBarry Smith     else {
57395373324SBarry Smith       /* assemble the entire matrix onto first processor. */
57495373324SBarry Smith       Mat     A;
5752ee70a88SLois Curfman McInnes       Mat_AIJ *Aloc;
5762eb8c8abSBarry Smith       int     M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
57795373324SBarry Smith       Scalar  *a;
5782ee70a88SLois Curfman McInnes 
57995373324SBarry Smith       if (!mytid) {
58095373324SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,0,0,0,&A);
58195373324SBarry Smith       }
58295373324SBarry Smith       else {
58395373324SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,0,0,0,&A);
58495373324SBarry Smith       }
58595373324SBarry Smith       CHKERR(ierr);
58695373324SBarry Smith 
58795373324SBarry Smith       /* copy over the A part */
5882ee70a88SLois Curfman McInnes       Aloc = (Mat_AIJ*) aij->A->data;
5892ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
59095373324SBarry Smith       row = aij->rstart;
59195373324SBarry Smith       for ( i=0; i<ai[m]; i++ ) {aj[i] += aij->cstart - 1;}
59295373324SBarry Smith       for ( i=0; i<m; i++ ) {
59307a0e7adSLois Curfman McInnes         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERTVALUES);
59495373324SBarry Smith         CHKERR(ierr);
59595373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
59695373324SBarry Smith       }
5972ee70a88SLois Curfman McInnes       aj = Aloc->j;
59895373324SBarry Smith       for ( i=0; i<ai[m]; i++ ) {aj[i] -= aij->cstart - 1;}
59995373324SBarry Smith 
60095373324SBarry Smith       /* copy over the B part */
6012ee70a88SLois Curfman McInnes       Aloc = (Mat_AIJ*) aij->B->data;
6022ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
60395373324SBarry Smith       row = aij->rstart;
60495373324SBarry Smith       ct = cols = (int *) MALLOC( (ai[m]+1)*sizeof(int) ); CHKPTR(cols);
60595373324SBarry Smith       for ( i=0; i<ai[m]; i++ ) {cols[i] = aij->garray[aj[i]-1];}
60695373324SBarry Smith       for ( i=0; i<m; i++ ) {
60707a0e7adSLois Curfman McInnes         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERTVALUES);
60895373324SBarry Smith         CHKERR(ierr);
60995373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
61095373324SBarry Smith       }
61195373324SBarry Smith       FREE(ct);
612ee50ffe9SBarry Smith       ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERR(ierr);
613ee50ffe9SBarry Smith       ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERR(ierr);
61495373324SBarry Smith       if (!mytid) {
61595373324SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERR(ierr);
61695373324SBarry Smith       }
61795373324SBarry Smith       ierr = MatDestroy(A); CHKERR(ierr);
61895373324SBarry Smith     }
61995373324SBarry Smith   }
6201eb62cbbSBarry Smith   return 0;
6211eb62cbbSBarry Smith }
6221eb62cbbSBarry Smith 
623d3c9fed9SLois Curfman McInnes extern int MatMarkDiag_AIJ(Mat_AIJ  *);
6241eb62cbbSBarry Smith /*
6251eb62cbbSBarry Smith     This has to provide several versions.
6261eb62cbbSBarry Smith 
6271eb62cbbSBarry Smith      1) per sequential
6281eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
6291eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
630d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
6311eb62cbbSBarry Smith */
632fc76ce87SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
633fc76ce87SLois Curfman McInnes                            double shift,int its,Vec xx)
6348a729477SBarry Smith {
63544a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
636d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
63744a69424SLois Curfman McInnes   Mat_AIJ    *A = (Mat_AIJ *) AA->data, *B = (Mat_AIJ *)BB->data;
6386abc6512SBarry Smith   Scalar     zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts;
6396abc6512SBarry Smith   int        ierr,*idx, *diag;
6406abc6512SBarry Smith   int        n = mat->n, m = mat->m, i;
641da3a660dSBarry Smith   Vec        tt;
6428a729477SBarry Smith 
64344a69424SLois Curfman McInnes   if (!mat->assembled) SETERR(1,"MatRelax_MPIAIJ: must assemble matrix first");
6441eb62cbbSBarry Smith 
645d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
646d6dfbf8fSBarry Smith   xs = x -1; /* shift by one for index start of 1 */
647da3a660dSBarry Smith   ls--;
648d3c9fed9SLois Curfman McInnes   if (!A->diag) {if ((ierr = MatMarkDiag_AIJ(A))) return ierr;}
649d6dfbf8fSBarry Smith   diag = A->diag;
650acb40c82SBarry Smith   if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) {
651acb40c82SBarry Smith     SETERR(1,"That option not yet support for parallel AIJ matrices");
652acb40c82SBarry Smith   }
653da3a660dSBarry Smith   if (flag & SOR_EISENSTAT) {
654da3a660dSBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
655da3a660dSBarry Smith     U is upper triangular, E is diagonal; This routine applies
656da3a660dSBarry Smith 
657da3a660dSBarry Smith             (L + E)^{-1} A (U + E)^{-1}
658da3a660dSBarry Smith 
659da3a660dSBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
660da3a660dSBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
661da3a660dSBarry Smith     is the relaxation factor.
662da3a660dSBarry Smith     */
663da3a660dSBarry Smith     ierr = VecCreate(xx,&tt); CHKERR(ierr);
664da3a660dSBarry Smith     VecGetArray(tt,&t);
665da3a660dSBarry Smith     scale = (2.0/omega) - 1.0;
666da3a660dSBarry Smith     /*  x = (E + U)^{-1} b */
667da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
668*9d00d63dSBarry Smith     ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
669da3a660dSBarry Smith                               mat->Mvctx); CHKERR(ierr);
670da3a660dSBarry Smith     for ( i=m-1; i>-1; i-- ) {
671da3a660dSBarry Smith       n    = A->i[i+1] - diag[i] - 1;
672da3a660dSBarry Smith       idx  = A->j + diag[i];
673da3a660dSBarry Smith       v    = A->a + diag[i];
674da3a660dSBarry Smith       sum  = b[i];
675da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
676da3a660dSBarry Smith       d    = shift + A->a[diag[i]-1];
677da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
678da3a660dSBarry Smith       idx  = B->j + B->i[i] - 1;
679da3a660dSBarry Smith       v    = B->a + B->i[i] - 1;
680da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
681da3a660dSBarry Smith       x[i] = omega*(sum/d);
682da3a660dSBarry Smith     }
683*9d00d63dSBarry Smith     ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
684da3a660dSBarry Smith                             mat->Mvctx); CHKERR(ierr);
685da3a660dSBarry Smith 
686da3a660dSBarry Smith     /*  t = b - (2*E - D)x */
687da3a660dSBarry Smith     v = A->a;
688da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ - 1])*x[i]; }
689da3a660dSBarry Smith 
690da3a660dSBarry Smith     /*  t = (E + L)^{-1}t */
691da3a660dSBarry Smith     ts = t - 1; /* shifted by one for index start of a or mat->j*/
692da3a660dSBarry Smith     diag = A->diag;
693da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
694*9d00d63dSBarry Smith     ierr = VecPipelineBegin(tt,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
695da3a660dSBarry Smith                                                  mat->Mvctx); CHKERR(ierr);
696da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
697da3a660dSBarry Smith       n    = diag[i] - A->i[i];
698da3a660dSBarry Smith       idx  = A->j + A->i[i] - 1;
699da3a660dSBarry Smith       v    = A->a + A->i[i] - 1;
700da3a660dSBarry Smith       sum  = t[i];
701da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
702da3a660dSBarry Smith       d    = shift + A->a[diag[i]-1];
703da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
704da3a660dSBarry Smith       idx  = B->j + B->i[i] - 1;
705da3a660dSBarry Smith       v    = B->a + B->i[i] - 1;
706da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
707da3a660dSBarry Smith       t[i] = omega*(sum/d);
708da3a660dSBarry Smith     }
709*9d00d63dSBarry Smith     ierr = VecPipelineEnd(tt,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
710da3a660dSBarry Smith                                                     mat->Mvctx); CHKERR(ierr);
711da3a660dSBarry Smith     /*  x = x + t */
712da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
713da3a660dSBarry Smith     VecDestroy(tt);
714da3a660dSBarry Smith     return 0;
715da3a660dSBarry Smith   }
716da3a660dSBarry Smith 
7171eb62cbbSBarry Smith 
718d6dfbf8fSBarry Smith   if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){
719da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
720da3a660dSBarry Smith       VecSet(&zero,mat->lvec); VecSet(&zero,xx);
721da3a660dSBarry Smith     }
722da3a660dSBarry Smith     else {
723*9d00d63dSBarry Smith       ierr=VecScatterBegin(xx,0,mat->lvec,0,INSERTVALUES,SCATTERUP,mat->Mvctx);
724d6dfbf8fSBarry Smith       CHKERR(ierr);
725*9d00d63dSBarry Smith       ierr = VecScatterEnd(xx,0,mat->lvec,0,INSERTVALUES,SCATTERUP,mat->Mvctx);
726d6dfbf8fSBarry Smith       CHKERR(ierr);
727da3a660dSBarry Smith     }
728d6dfbf8fSBarry Smith     while (its--) {
729d6dfbf8fSBarry Smith       /* go down through the rows */
730*9d00d63dSBarry Smith       ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
731d6dfbf8fSBarry Smith                               mat->Mvctx); CHKERR(ierr);
732d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
733d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
734d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
735d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
736d6dfbf8fSBarry Smith         sum  = b[i];
737d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
738d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
739d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
740d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
741d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
742d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
743d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
744d6dfbf8fSBarry Smith       }
745*9d00d63dSBarry Smith       ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
746d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
747d6dfbf8fSBarry Smith       /* come up through the rows */
748*9d00d63dSBarry Smith       ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
749d6dfbf8fSBarry Smith                               mat->Mvctx); CHKERR(ierr);
750d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
751d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
752d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
753d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
754d6dfbf8fSBarry Smith         sum  = b[i];
755d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
756d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
757d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
758d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
759d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
760d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
761d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
762d6dfbf8fSBarry Smith       }
763*9d00d63dSBarry Smith       ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
764d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
765d6dfbf8fSBarry Smith     }
766d6dfbf8fSBarry Smith   }
767d6dfbf8fSBarry Smith   else if (flag & SOR_FORWARD_SWEEP){
768da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
769da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
770*9d00d63dSBarry Smith       ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
771da3a660dSBarry Smith                               mat->Mvctx); CHKERR(ierr);
772da3a660dSBarry Smith       for ( i=0; i<m; i++ ) {
773da3a660dSBarry Smith         n    = diag[i] - A->i[i];
774da3a660dSBarry Smith         idx  = A->j + A->i[i] - 1;
775da3a660dSBarry Smith         v    = A->a + A->i[i] - 1;
776da3a660dSBarry Smith         sum  = b[i];
777da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
778da3a660dSBarry Smith         d    = shift + A->a[diag[i]-1];
779da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
780da3a660dSBarry Smith         idx  = B->j + B->i[i] - 1;
781da3a660dSBarry Smith         v    = B->a + B->i[i] - 1;
782da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
783da3a660dSBarry Smith         x[i] = omega*(sum/d);
784da3a660dSBarry Smith       }
785*9d00d63dSBarry Smith       ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
786da3a660dSBarry Smith                             mat->Mvctx); CHKERR(ierr);
787da3a660dSBarry Smith       its--;
788da3a660dSBarry Smith     }
789d6dfbf8fSBarry Smith     while (its--) {
790*9d00d63dSBarry Smith       ierr=VecScatterBegin(xx,0,mat->lvec,0,INSERTVALUES,SCATTERUP,mat->Mvctx);
791d6dfbf8fSBarry Smith       CHKERR(ierr);
792*9d00d63dSBarry Smith       ierr = VecScatterEnd(xx,0,mat->lvec,0,INSERTVALUES,SCATTERUP,mat->Mvctx);
793d6dfbf8fSBarry Smith       CHKERR(ierr);
794*9d00d63dSBarry Smith       ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
795d6dfbf8fSBarry Smith                               mat->Mvctx); CHKERR(ierr);
796d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
797d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
798d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
799d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
800d6dfbf8fSBarry Smith         sum  = b[i];
801d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
802d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
803d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
804d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
805d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
806d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
807d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
808d6dfbf8fSBarry Smith       }
809*9d00d63dSBarry Smith       ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEDOWN,
810d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
811d6dfbf8fSBarry Smith     }
812d6dfbf8fSBarry Smith   }
813d6dfbf8fSBarry Smith   else if (flag & SOR_BACKWARD_SWEEP){
814da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
815da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
816*9d00d63dSBarry Smith       ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
817da3a660dSBarry Smith                               mat->Mvctx); CHKERR(ierr);
818da3a660dSBarry Smith       for ( i=m-1; i>-1; i-- ) {
819da3a660dSBarry Smith         n    = A->i[i+1] - diag[i] - 1;
820da3a660dSBarry Smith         idx  = A->j + diag[i];
821da3a660dSBarry Smith         v    = A->a + diag[i];
822da3a660dSBarry Smith         sum  = b[i];
823da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
824da3a660dSBarry Smith         d    = shift + A->a[diag[i]-1];
825da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
826da3a660dSBarry Smith         idx  = B->j + B->i[i] - 1;
827da3a660dSBarry Smith         v    = B->a + B->i[i] - 1;
828da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
829da3a660dSBarry Smith         x[i] = omega*(sum/d);
830da3a660dSBarry Smith       }
831*9d00d63dSBarry Smith       ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
832da3a660dSBarry Smith                             mat->Mvctx); CHKERR(ierr);
833da3a660dSBarry Smith       its--;
834da3a660dSBarry Smith     }
835d6dfbf8fSBarry Smith     while (its--) {
836*9d00d63dSBarry Smith       ierr = VecScatterBegin(xx,0,mat->lvec,0,INSERTVALUES,SCATTERDOWN,
837d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
838*9d00d63dSBarry Smith       ierr = VecScatterEnd(xx,0,mat->lvec,0,INSERTVALUES,SCATTERDOWN,
839d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
840*9d00d63dSBarry Smith       ierr = VecPipelineBegin(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
841d6dfbf8fSBarry Smith                               mat->Mvctx); CHKERR(ierr);
842d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
843d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
844d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
845d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
846d6dfbf8fSBarry Smith         sum  = b[i];
847d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
848d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
849d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
850d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
851d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
852d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
853d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
854d6dfbf8fSBarry Smith       }
855*9d00d63dSBarry Smith       ierr = VecPipelineEnd(xx,0,mat->lvec,0,INSERTVALUES,PIPELINEUP,
856d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
857d6dfbf8fSBarry Smith     }
858d6dfbf8fSBarry Smith   }
859d6dfbf8fSBarry Smith   else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
860da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
861da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
862da3a660dSBarry Smith     }
863*9d00d63dSBarry Smith     ierr=VecScatterBegin(xx,0,mat->lvec,0,INSERTVALUES,SCATTERALL,mat->Mvctx);
864d6dfbf8fSBarry Smith     CHKERR(ierr);
865*9d00d63dSBarry Smith     ierr = VecScatterEnd(xx,0,mat->lvec,0,INSERTVALUES,SCATTERALL,mat->Mvctx);
866d6dfbf8fSBarry Smith     CHKERR(ierr);
867d6dfbf8fSBarry Smith     while (its--) {
868d6dfbf8fSBarry Smith       /* go down through the rows */
869d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
870d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
871d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
872d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
873d6dfbf8fSBarry Smith         sum  = b[i];
874d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
875d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
876d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
877d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
878d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
879d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
880d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
881d6dfbf8fSBarry Smith       }
882d6dfbf8fSBarry Smith       /* come up through the rows */
883d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
884d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
885d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
886d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
887d6dfbf8fSBarry Smith         sum  = b[i];
888d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
889d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
890d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
891d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
892d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
893d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
894d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
895d6dfbf8fSBarry Smith       }
896d6dfbf8fSBarry Smith     }
897d6dfbf8fSBarry Smith   }
898d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_FORWARD_SWEEP){
899da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
900da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
901da3a660dSBarry Smith     }
902*9d00d63dSBarry Smith     ierr=VecScatterBegin(xx,0,mat->lvec,0,INSERTVALUES,SCATTERALL,mat->Mvctx);
903d6dfbf8fSBarry Smith     CHKERR(ierr);
904*9d00d63dSBarry Smith     ierr = VecScatterEnd(xx,0,mat->lvec,0,INSERTVALUES,SCATTERALL,mat->Mvctx);
905d6dfbf8fSBarry Smith     CHKERR(ierr);
906d6dfbf8fSBarry Smith     while (its--) {
907d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
908d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
909d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
910d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
911d6dfbf8fSBarry Smith         sum  = b[i];
912d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
913d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
914d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
915d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
916d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
917d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
918d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
919d6dfbf8fSBarry Smith       }
920d6dfbf8fSBarry Smith     }
921d6dfbf8fSBarry Smith   }
922d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_BACKWARD_SWEEP){
923da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
924da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
925da3a660dSBarry Smith     }
926*9d00d63dSBarry Smith     ierr = VecScatterBegin(xx,0,mat->lvec,0,INSERTVALUES,SCATTERALL,
927d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
928*9d00d63dSBarry Smith     ierr = VecScatterEnd(xx,0,mat->lvec,0,INSERTVALUES,SCATTERALL,
929d6dfbf8fSBarry Smith                             mat->Mvctx); CHKERR(ierr);
930d6dfbf8fSBarry Smith     while (its--) {
931d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
932d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
933d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
934d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
935d6dfbf8fSBarry Smith         sum  = b[i];
936d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
937d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
938d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
939d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
940d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
941d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
942d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
943d6dfbf8fSBarry Smith       }
944d6dfbf8fSBarry Smith     }
945d6dfbf8fSBarry Smith   }
9468a729477SBarry Smith   return 0;
9478a729477SBarry Smith }
948a66be287SLois Curfman McInnes 
949fc76ce87SLois Curfman McInnes static int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,int *nz,
950fc76ce87SLois Curfman McInnes                              int *nzalloc,int *mem)
951a66be287SLois Curfman McInnes {
952a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
953a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
954a66be287SLois Curfman McInnes   int        ierr, isend[3], irecv[3], nzA, nzallocA, memA;
955a66be287SLois Curfman McInnes 
956a66be287SLois Curfman McInnes   ierr = MatGetInfo(A,MAT_LOCAL,&nzA,&nzallocA,&memA); CHKERR(ierr);
957a66be287SLois Curfman McInnes   ierr = MatGetInfo(B,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERR(ierr);
958a66be287SLois Curfman McInnes   isend[0] += nzA; isend[1] += nzallocA; isend[2] += memA;
959a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
960a66be287SLois Curfman McInnes     *nz = isend[0]; *nzalloc = isend[1]; *mem = isend[2];
961a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
962a66be287SLois Curfman McInnes     MPI_Allreduce((void *) isend,(void *) irecv,3,MPI_INT,MPI_MAX,matin->comm);
963a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
964a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
965a66be287SLois Curfman McInnes     MPI_Allreduce((void *) isend,(void *) irecv,3,MPI_INT,MPI_SUM,matin->comm);
966a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
967a66be287SLois Curfman McInnes   }
968a66be287SLois Curfman McInnes   return 0;
969a66be287SLois Curfman McInnes }
970a66be287SLois Curfman McInnes 
971fc76ce87SLois Curfman McInnes static int MatSetOption_MPIAIJ(Mat aijin,MatOption op)
972c74985f6SBarry Smith {
97344a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
974c74985f6SBarry Smith 
975c74985f6SBarry Smith   if      (op == NO_NEW_NONZERO_LOCATIONS)  {
976c74985f6SBarry Smith     MatSetOption(aij->A,op);
977c74985f6SBarry Smith     MatSetOption(aij->B,op);
978c74985f6SBarry Smith   }
979c74985f6SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) {
980c74985f6SBarry Smith     MatSetOption(aij->A,op);
981c74985f6SBarry Smith     MatSetOption(aij->B,op);
982c74985f6SBarry Smith   }
983c74985f6SBarry Smith   else if (op == COLUMN_ORIENTED) SETERR(1,"Column oriented not supported");
984c74985f6SBarry Smith   return 0;
985c74985f6SBarry Smith }
986c74985f6SBarry Smith 
987d1710a03SLois Curfman McInnes static int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
988c74985f6SBarry Smith {
98944a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
990c74985f6SBarry Smith   *m = mat->M; *n = mat->N;
991c74985f6SBarry Smith   return 0;
992c74985f6SBarry Smith }
993c74985f6SBarry Smith 
994d1710a03SLois Curfman McInnes static int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
995c74985f6SBarry Smith {
99644a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
997c74985f6SBarry Smith   *m = mat->m; *n = mat->n;
998c74985f6SBarry Smith   return 0;
999c74985f6SBarry Smith }
1000c74985f6SBarry Smith 
1001d1710a03SLois Curfman McInnes static int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
1002c74985f6SBarry Smith {
100344a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1004c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
1005c74985f6SBarry Smith   return 0;
1006c74985f6SBarry Smith }
1007c74985f6SBarry Smith 
100839e00950SLois Curfman McInnes static int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
100939e00950SLois Curfman McInnes {
1010154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1011154123eaSLois Curfman McInnes   Scalar     *vworkA, *vworkB, **pvA, **pvB;
1012154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
1013154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
101439e00950SLois Curfman McInnes 
1015154123eaSLois Curfman McInnes   if (!mat->assembled)
101639e00950SLois Curfman McInnes     SETERR(1,"MatGetRow_MPIAIJ: Must assemble matrix first.");
101739e00950SLois Curfman McInnes   if (row < rstart || row >= rend)
1018abc0e9e4SLois Curfman McInnes     SETERR(1,"MatGetRow_MPIAIJ: Currently you can get only local rows.")
1019abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
102039e00950SLois Curfman McInnes 
1021154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1022154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1023154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
1024154123eaSLois Curfman McInnes   ierr = MatGetRow(mat->A,lrow,&nzA,pcA,pvA); CHKERR(ierr);
1025154123eaSLois Curfman McInnes   ierr = MatGetRow(mat->B,lrow,&nzB,pcB,pvB); CHKERR(ierr);
1026154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1027154123eaSLois Curfman McInnes 
1028154123eaSLois Curfman McInnes   if (v  || idx) {
1029154123eaSLois Curfman McInnes     if (nztot) {
1030154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
1031154123eaSLois Curfman McInnes       int imark, imark2;
1032154123eaSLois Curfman McInnes       for (i=0; i<nzB; i++) cworkB[i] = mat->garray[cworkB[i]];
1033154123eaSLois Curfman McInnes       if (v) {
103439e00950SLois Curfman McInnes         *v = (Scalar *) MALLOC( (nztot)*sizeof(Scalar) ); CHKPTR(*v);
103539e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1036154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*v)[i] = vworkB[i];
1037154123eaSLois Curfman McInnes           else break;
1038154123eaSLois Curfman McInnes         }
1039154123eaSLois Curfman McInnes         imark = i;
1040154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*v)[imark+i] = vworkA[i];
1041154123eaSLois Curfman McInnes         imark2 = imark+nzA;
1042154123eaSLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*v)[imark2+i] = vworkB[i];
1043154123eaSLois Curfman McInnes       }
1044154123eaSLois Curfman McInnes       if (idx) {
1045154123eaSLois Curfman McInnes         *idx = (int *) MALLOC( (nztot)*sizeof(int) ); CHKPTR(*idx);
1046154123eaSLois Curfman McInnes         for (i=0; i<nzA; i++) cworkA[i] += cstart;
1047154123eaSLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1048154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*idx)[i] = cworkB[i];
1049154123eaSLois Curfman McInnes           else break;
1050154123eaSLois Curfman McInnes         }
1051154123eaSLois Curfman McInnes         imark = i;
1052154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*idx)[imark+i] = cworkA[i];
1053154123eaSLois Curfman McInnes         imark2 = imark+nzA;
1054154123eaSLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*idx)[imark2+i] = cworkB[i];
105539e00950SLois Curfman McInnes       }
105639e00950SLois Curfman McInnes     }
105739e00950SLois Curfman McInnes     else {*idx = 0; *v=0;}
1058154123eaSLois Curfman McInnes   }
105939e00950SLois Curfman McInnes   *nz = nztot;
1060154123eaSLois Curfman McInnes   ierr = MatRestoreRow(mat->A,lrow,&nzA,pcA,pvA); CHKERR(ierr);
1061154123eaSLois Curfman McInnes   ierr = MatRestoreRow(mat->B,lrow,&nzB,pcB,pvB); CHKERR(ierr);
106239e00950SLois Curfman McInnes   return 0;
106339e00950SLois Curfman McInnes }
106439e00950SLois Curfman McInnes 
106539e00950SLois Curfman McInnes static int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
106639e00950SLois Curfman McInnes {
1067154123eaSLois Curfman McInnes   if (idx) FREE(*idx);
1068154123eaSLois Curfman McInnes   if (v) FREE(*v);
106939e00950SLois Curfman McInnes   return 0;
107039e00950SLois Curfman McInnes }
107139e00950SLois Curfman McInnes 
1072fc76ce87SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MatType,Mat *);
107307a0e7adSLois Curfman McInnes static int MatCopy_MPIAIJ_Private(Mat,Mat *);
1074d6dfbf8fSBarry Smith 
10758a729477SBarry Smith /* -------------------------------------------------------------------*/
10762ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
107739e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
107844a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
107944a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
10801eb62cbbSBarry Smith        0,0,0,0,
10818a729477SBarry Smith        0,0,
108244a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
10838a729477SBarry Smith        0,
1084a66be287SLois Curfman McInnes        MatGetInfo_MPIAIJ,0,
1085d1710a03SLois Curfman McInnes        MatGetDiagonal_MPIAIJ,0,0,
1086ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
10871eb62cbbSBarry Smith        0,
10882ee70a88SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,0,
1089c74985f6SBarry Smith        0,0,0,0,
1090d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
109177915d1cSLois Curfman McInnes        0,0,
109207a0e7adSLois Curfman McInnes        0,0,MatConvert_MPIAIJ,0,0,MatCopy_MPIAIJ_Private};
10938a729477SBarry Smith 
10948a729477SBarry Smith /*@
1095ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
1096ff756334SLois Curfman McInnes    (the default parallel PETSc format).
10978a729477SBarry Smith 
10988a729477SBarry Smith    Input Parameters:
10991eb62cbbSBarry Smith .  comm - MPI communicator
11007d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
11017d3e4905SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated
11027d3e4905SLois Curfman McInnes            if N is given)
11037d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
11047d3e4905SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated
11057d3e4905SLois Curfman McInnes            if n is given)
1106ff756334SLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of matrix
1107ff756334SLois Curfman McInnes            (same for all local rows)
11081eb62cbbSBarry Smith .  d_nzz - number of nonzeros per row in diagonal portion of matrix or null
1109ff756334SLois Curfman McInnes            (possibly different for each row).  You must leave room for the
1110ff756334SLois Curfman McInnes            diagonal entry even if it is zero.
1111ff756334SLois Curfman McInnes .  o_nz - number of nonzeros per row in off-diagonal portion of matrix
1112ff756334SLois Curfman McInnes            (same for all local rows)
11131eb62cbbSBarry Smith .  o_nzz - number of nonzeros per row in off-diagonal portion of matrix
1114ff756334SLois Curfman McInnes            or null (possibly different for each row).
11158a729477SBarry Smith 
1116ff756334SLois Curfman McInnes    Output Parameter:
11178a729477SBarry Smith .  newmat - the matrix
11188a729477SBarry Smith 
1119ff756334SLois Curfman McInnes    Notes:
1120ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1121ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
1122ff756334SLois Curfman McInnes    storage.  That is, the stored row and column indices begin at
1123ff756334SLois Curfman McInnes    one, not zero.
1124ff756334SLois Curfman McInnes 
1125ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1126ff756334SLois Curfman McInnes    (possibly both).
1127ff756334SLois Curfman McInnes 
1128ff756334SLois Curfman McInnes    The user can set d_nz, d_nnz, o_nz, and o_nnz to zero for PETSc to
1129ff756334SLois Curfman McInnes    control dynamic memory allocation.
1130ff756334SLois Curfman McInnes 
1131ff756334SLois Curfman McInnes .keywords: Mat, matrix, aij, compressed row, sparse, parallel
1132ff756334SLois Curfman McInnes 
1133ff756334SLois Curfman McInnes .seealso: MatCreateSequentialAIJ(), MatSetValues()
11348a729477SBarry Smith @*/
11351eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,
11361eb62cbbSBarry Smith                  int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat)
11378a729477SBarry Smith {
11388a729477SBarry Smith   Mat          mat;
113944a69424SLois Curfman McInnes   Mat_MPIAIJ   *aij;
11406abc6512SBarry Smith   int          ierr, i,sum[2],work[2];
11418a729477SBarry Smith   *newmat         = 0;
114244a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm);
1143a5a9c739SBarry Smith   PLogObjectCreate(mat);
114444a69424SLois Curfman McInnes   mat->data       = (void *) (aij = NEW(Mat_MPIAIJ)); CHKPTR(aij);
11458a729477SBarry Smith   mat->ops        = &MatOps;
114644a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
114744a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
11488a729477SBarry Smith   mat->factor     = 0;
1149d6dfbf8fSBarry Smith 
1150d6dfbf8fSBarry Smith   mat->comm       = comm;
115107a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
11521eb62cbbSBarry Smith   MPI_Comm_rank(comm,&aij->mytid);
11531eb62cbbSBarry Smith   MPI_Comm_size(comm,&aij->numtids);
11541eb62cbbSBarry Smith 
11551eb62cbbSBarry Smith   if (M == -1 || N == -1) {
11561eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1157d6dfbf8fSBarry Smith     MPI_Allreduce((void *) work,(void *) sum,2,MPI_INT,MPI_SUM,comm );
11581eb62cbbSBarry Smith     if (M == -1) M = sum[0];
11591eb62cbbSBarry Smith     if (N == -1) N = sum[1];
11601eb62cbbSBarry Smith   }
11611eb62cbbSBarry Smith   if (m == -1) {m = M/aij->numtids + ((M % aij->numtids) > aij->mytid);}
11621eb62cbbSBarry Smith   if (n == -1) {n = N/aij->numtids + ((N % aij->numtids) > aij->mytid);}
11638a729477SBarry Smith   aij->m       = m;
11648a729477SBarry Smith   aij->n       = n;
11651eb62cbbSBarry Smith   aij->N       = N;
11661eb62cbbSBarry Smith   aij->M       = M;
11671eb62cbbSBarry Smith 
11681eb62cbbSBarry Smith   /* build local table of row and column ownerships */
11691eb62cbbSBarry Smith   aij->rowners = (int *) MALLOC(2*(aij->numtids+2)*sizeof(int));
11701eb62cbbSBarry Smith   CHKPTR(aij->rowners);
11711eb62cbbSBarry Smith   aij->cowners = aij->rowners + aij->numtids + 1;
11721eb62cbbSBarry Smith   MPI_Allgather(&m,1,MPI_INT,aij->rowners+1,1,MPI_INT,comm);
11731eb62cbbSBarry Smith   aij->rowners[0] = 0;
11741eb62cbbSBarry Smith   for ( i=2; i<=aij->numtids; i++ ) {
11751eb62cbbSBarry Smith     aij->rowners[i] += aij->rowners[i-1];
11768a729477SBarry Smith   }
11771eb62cbbSBarry Smith   aij->rstart = aij->rowners[aij->mytid];
11781eb62cbbSBarry Smith   aij->rend   = aij->rowners[aij->mytid+1];
11791eb62cbbSBarry Smith   MPI_Allgather(&n,1,MPI_INT,aij->cowners+1,1,MPI_INT,comm);
11801eb62cbbSBarry Smith   aij->cowners[0] = 0;
11811eb62cbbSBarry Smith   for ( i=2; i<=aij->numtids; i++ ) {
11821eb62cbbSBarry Smith     aij->cowners[i] += aij->cowners[i-1];
11838a729477SBarry Smith   }
11841eb62cbbSBarry Smith   aij->cstart = aij->cowners[aij->mytid];
11851eb62cbbSBarry Smith   aij->cend   = aij->cowners[aij->mytid+1];
11868a729477SBarry Smith 
11878a729477SBarry Smith 
11886b5873e3SBarry Smith   ierr = MatCreateSequentialAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&aij->A);
11896b5873e3SBarry Smith   CHKERR(ierr);
1190a5a9c739SBarry Smith   PLogObjectParent(mat,aij->A);
11916b5873e3SBarry Smith   ierr = MatCreateSequentialAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&aij->B);
11926b5873e3SBarry Smith   CHKERR(ierr);
1193a5a9c739SBarry Smith   PLogObjectParent(mat,aij->B);
11948a729477SBarry Smith 
11951eb62cbbSBarry Smith   /* build cache for off array entries formed */
11961eb62cbbSBarry Smith   aij->stash.nmax = CHUNCKSIZE; /* completely arbratray number */
11971eb62cbbSBarry Smith   aij->stash.n    = 0;
11981eb62cbbSBarry Smith   aij->stash.array = (Scalar *) MALLOC( aij->stash.nmax*(2*sizeof(int) +
11991eb62cbbSBarry Smith                             sizeof(Scalar))); CHKPTR(aij->stash.array);
12001eb62cbbSBarry Smith   aij->stash.idx = (int *) (aij->stash.array + aij->stash.nmax);
12011eb62cbbSBarry Smith   aij->stash.idy = (int *) (aij->stash.idx + aij->stash.nmax);
12029e25ed09SBarry Smith   aij->colmap    = 0;
12039e25ed09SBarry Smith   aij->garray    = 0;
12048a729477SBarry Smith 
12051eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
12061eb62cbbSBarry Smith   aij->lvec      = 0;
12071eb62cbbSBarry Smith   aij->Mvctx     = 0;
1208d6dfbf8fSBarry Smith   aij->assembled = 0;
12098a729477SBarry Smith 
1210d6dfbf8fSBarry Smith   *newmat = mat;
1211d6dfbf8fSBarry Smith   return 0;
1212d6dfbf8fSBarry Smith }
1213c74985f6SBarry Smith 
121407a0e7adSLois Curfman McInnes static int MatCopy_MPIAIJ_Private(Mat matin,Mat *newmat)
1215d6dfbf8fSBarry Smith {
1216d6dfbf8fSBarry Smith   Mat        mat;
121744a69424SLois Curfman McInnes   Mat_MPIAIJ *aij,*oldmat = (Mat_MPIAIJ *) matin->data;
12182ee70a88SLois Curfman McInnes   int        ierr, len;
1219d6dfbf8fSBarry Smith   *newmat      = 0;
1220d6dfbf8fSBarry Smith 
1221d6dfbf8fSBarry Smith   if (!oldmat->assembled) SETERR(1,"Cannot copy unassembled matrix");
122244a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm);
1223a5a9c739SBarry Smith   PLogObjectCreate(mat);
122444a69424SLois Curfman McInnes   mat->data       = (void *) (aij = NEW(Mat_MPIAIJ)); CHKPTR(aij);
1225d6dfbf8fSBarry Smith   mat->ops        = &MatOps;
122644a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
122744a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
1228d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1229d6dfbf8fSBarry Smith 
1230d6dfbf8fSBarry Smith   aij->m          = oldmat->m;
1231d6dfbf8fSBarry Smith   aij->n          = oldmat->n;
1232d6dfbf8fSBarry Smith   aij->M          = oldmat->M;
1233d6dfbf8fSBarry Smith   aij->N          = oldmat->N;
1234d6dfbf8fSBarry Smith 
1235d6dfbf8fSBarry Smith   aij->assembled  = 1;
1236d6dfbf8fSBarry Smith   aij->rstart     = oldmat->rstart;
1237d6dfbf8fSBarry Smith   aij->rend       = oldmat->rend;
1238d6dfbf8fSBarry Smith   aij->cstart     = oldmat->cstart;
1239d6dfbf8fSBarry Smith   aij->cend       = oldmat->cend;
1240d6dfbf8fSBarry Smith   aij->numtids    = oldmat->numtids;
1241d6dfbf8fSBarry Smith   aij->mytid      = oldmat->mytid;
124207a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
1243d6dfbf8fSBarry Smith 
1244d6dfbf8fSBarry Smith   aij->rowners        = (int *) MALLOC( (aij->numtids+1)*sizeof(int) );
1245d6dfbf8fSBarry Smith   CHKPTR(aij->rowners);
1246d6dfbf8fSBarry Smith   MEMCPY(aij->rowners,oldmat->rowners,(aij->numtids+1)*sizeof(int));
1247d6dfbf8fSBarry Smith   aij->stash.nmax     = 0;
1248d6dfbf8fSBarry Smith   aij->stash.n        = 0;
1249d6dfbf8fSBarry Smith   aij->stash.array    = 0;
12502ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
12512ee70a88SLois Curfman McInnes     aij->colmap      = (int *) MALLOC( (aij->N)*sizeof(int) );
12522ee70a88SLois Curfman McInnes     CHKPTR(aij->colmap);
12532ee70a88SLois Curfman McInnes     MEMCPY(aij->colmap,oldmat->colmap,(aij->N)*sizeof(int));
12542ee70a88SLois Curfman McInnes   } else aij->colmap = 0;
12552ee70a88SLois Curfman McInnes     if (oldmat->garray && (len = ((Mat_AIJ *) (oldmat->B->data))->n)) {
12562ee70a88SLois Curfman McInnes     aij->garray      = (int *) MALLOC(len*sizeof(int) ); CHKPTR(aij->garray);
12572ee70a88SLois Curfman McInnes     MEMCPY(aij->garray,oldmat->garray,len*sizeof(int));
12582ee70a88SLois Curfman McInnes   } else aij->garray = 0;
1259d6dfbf8fSBarry Smith   mat->comm           = matin->comm;
1260d6dfbf8fSBarry Smith 
1261d6dfbf8fSBarry Smith   ierr =  VecCreate(oldmat->lvec,&aij->lvec); CHKERR(ierr);
1262a5a9c739SBarry Smith   PLogObjectParent(mat,aij->lvec);
1263d6dfbf8fSBarry Smith   ierr =  VecScatterCtxCopy(oldmat->Mvctx,&aij->Mvctx); CHKERR(ierr);
1264a5a9c739SBarry Smith   PLogObjectParent(mat,aij->Mvctx);
126507a0e7adSLois Curfman McInnes   ierr =  MatConvert(oldmat->A,MATSAME,&aij->A); CHKERR(ierr);
1266a5a9c739SBarry Smith   PLogObjectParent(mat,aij->A);
126707a0e7adSLois Curfman McInnes   ierr =  MatConvert(oldmat->B,MATSAME,&aij->B); CHKERR(ierr);
1268a5a9c739SBarry Smith   PLogObjectParent(mat,aij->B);
12698a729477SBarry Smith   *newmat = mat;
12708a729477SBarry Smith   return 0;
12718a729477SBarry Smith }
1272