xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision fafbff5397be748e60f6fa9fc551a26f19f3b045)
1cb512458SBarry Smith #ifndef lint
2*fafbff53SBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.77 1995/09/11 16:52:56 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 
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;
1744a69424SLois Curfman McInnes   Mat_AIJ    *B = (Mat_AIJ*) aij->B->data;
189e25ed09SBarry Smith   int        n = B->n,i;
1978b31e54SBarry Smith   aij->colmap = (int *) PETSCMALLOC(aij->N*sizeof(int));CHKPTRQ(aij->colmap);
20464493b3SBarry Smith   PLogObjectMemory(mat,aij->N*sizeof(int));
2178b31e54SBarry Smith   PETSCMEMSET(aij->colmap,0,aij->N*sizeof(int));
22abc0e9e4SLois Curfman McInnes   for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1;
239e25ed09SBarry Smith   return 0;
249e25ed09SBarry Smith }
259e25ed09SBarry Smith 
262493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat);
272493cbb0SBarry Smith 
2851c98ccdSLois Curfman McInnes static int MatGetReordering_MPIAIJ(Mat mat,MatOrdering type,IS *rperm,IS *cperm)
29299609e3SLois Curfman McInnes {
30299609e3SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
31299609e3SLois Curfman McInnes   int ierr;
32299609e3SLois Curfman McInnes   if (aij->numtids == 1) {
3351c98ccdSLois Curfman McInnes     ierr = MatGetReordering(aij->A,type,rperm,cperm); CHKERRQ(ierr);
34299609e3SLois Curfman McInnes   } else
35bbb6d6a8SBarry Smith     SETERRQ(1,"MatGetReordering_MPIAIJ:  not yet 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;
431eb62cbbSBarry Smith   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
441eb62cbbSBarry Smith   int        cstart = aij->cstart, cend = aij->cend,row,col;
458a729477SBarry Smith 
4607a0e7adSLois Curfman McInnes   if (aij->insertmode != NOTSETVALUES && aij->insertmode != addv) {
47bbb6d6a8SBarry Smith     SETERRQ(1,"MatSetValues_MPIAIJ:You cannot mix inserts and adds");
488a729477SBarry Smith   }
491eb62cbbSBarry Smith   aij->insertmode = addv;
508a729477SBarry Smith   for ( i=0; i<m; i++ ) {
51bbb6d6a8SBarry Smith     if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIAIJ:Negative row");
52bbb6d6a8SBarry Smith     if (idxm[i] >= aij->M) SETERRQ(1,"MatSetValues_MPIAIJ:Row too large");
531eb62cbbSBarry Smith     if (idxm[i] >= rstart && idxm[i] < rend) {
541eb62cbbSBarry Smith       row = idxm[i] - rstart;
551eb62cbbSBarry Smith       for ( j=0; j<n; j++ ) {
56bbb6d6a8SBarry Smith         if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIAIJ:Negative column");
57bbb6d6a8SBarry Smith         if (idxn[j] >= aij->N) SETERRQ(1,"MatSetValues_MPIAIJ:Col too large");
581eb62cbbSBarry Smith         if (idxn[j] >= cstart && idxn[j] < cend){
591eb62cbbSBarry Smith           col = idxn[j] - cstart;
6078b31e54SBarry Smith           ierr = MatSetValues(aij->A,1,&row,1,&col,v+i*n+j,addv);CHKERRQ(ierr);
611eb62cbbSBarry Smith         }
621eb62cbbSBarry Smith         else {
63d6dfbf8fSBarry Smith           if (aij->assembled) {
64b7c46309SBarry Smith             if (!aij->colmap) {ierr = CreateColmap_Private(mat);CHKERRQ(ierr);}
659e25ed09SBarry Smith             col = aij->colmap[idxn[j]] - 1;
66b7c46309SBarry Smith             if (col < 0 && !((Mat_AIJ*)(aij->A->data))->nonew) {
672493cbb0SBarry Smith               ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
682493cbb0SBarry Smith               col =  idxn[j];
69d6dfbf8fSBarry Smith             }
709e25ed09SBarry Smith           }
719e25ed09SBarry Smith           else col = idxn[j];
7278b31e54SBarry Smith           ierr = MatSetValues(aij->B,1,&row,1,&col,v+i*n+j,addv);CHKERRQ(ierr);
731eb62cbbSBarry Smith         }
741eb62cbbSBarry Smith       }
751eb62cbbSBarry Smith     }
761eb62cbbSBarry Smith     else {
77dbd7a890SLois Curfman McInnes       ierr = StashValues_Private(&aij->stash,idxm[i],n,idxn,v+i*n,addv);
7878b31e54SBarry Smith       CHKERRQ(ierr);
791eb62cbbSBarry Smith     }
808a729477SBarry Smith   }
818a729477SBarry Smith   return 0;
828a729477SBarry Smith }
838a729477SBarry Smith 
848a729477SBarry Smith /*
851eb62cbbSBarry Smith     the assembly code is a lot like the code for vectors, we should
861eb62cbbSBarry Smith     sometime derive a single assembly code that can be used for
871eb62cbbSBarry Smith     either case.
888a729477SBarry Smith */
898a729477SBarry Smith 
90fc76ce87SLois Curfman McInnes static int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode)
918a729477SBarry Smith {
9244a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
93d6dfbf8fSBarry Smith   MPI_Comm    comm = mat->comm;
946abc6512SBarry Smith   int         numtids = aij->numtids, *owners = aij->rowners;
951eb62cbbSBarry Smith   int         mytid = aij->mytid;
961eb62cbbSBarry Smith   MPI_Request *send_waits,*recv_waits;
976abc6512SBarry Smith   int         *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
985392566eSBarry Smith   int         tag = mat->tag, *owner,*starts,count,ierr;
991eb62cbbSBarry Smith   InsertMode  addv;
1001eb62cbbSBarry Smith   Scalar      *rvalues,*svalues;
1011eb62cbbSBarry Smith 
1021eb62cbbSBarry Smith   /* make sure all processors are either in INSERTMODE or ADDMODE */
10328988994SBarry Smith   MPI_Allreduce((void *) &aij->insertmode,(void *) &addv,1,MPI_INT,
1041eb62cbbSBarry Smith                 MPI_BOR,comm);
1059d00d63dSBarry Smith   if (addv == (ADDVALUES|INSERTVALUES)) {
106bbb6d6a8SBarry Smith     SETERRQ(1,"MatAssemblyBegin_MPIAIJ:Some processors inserted others added");
1071eb62cbbSBarry Smith   }
1081eb62cbbSBarry Smith   aij->insertmode = addv; /* in case this processor had no cache */
1091eb62cbbSBarry Smith 
1101eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
11178b31e54SBarry Smith   nprocs = (int *) PETSCMALLOC( 2*numtids*sizeof(int) ); CHKPTRQ(nprocs);
11278b31e54SBarry Smith   PETSCMEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids;
11378b31e54SBarry Smith   owner = (int *) PETSCMALLOC( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
1141eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
1151eb62cbbSBarry Smith     idx = aij->stash.idx[i];
1161eb62cbbSBarry Smith     for ( j=0; j<numtids; j++ ) {
1171eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
1181eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
1198a729477SBarry Smith       }
1208a729477SBarry Smith     }
1218a729477SBarry Smith   }
1221eb62cbbSBarry Smith   nsends = 0;  for ( i=0; i<numtids; i++ ) { nsends += procs[i];}
1231eb62cbbSBarry Smith 
1241eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
12578b31e54SBarry Smith   work = (int *) PETSCMALLOC( numtids*sizeof(int) ); CHKPTRQ(work);
1261eb62cbbSBarry Smith   MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm);
1271eb62cbbSBarry Smith   nreceives = work[mytid];
1281eb62cbbSBarry Smith   MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm);
1291eb62cbbSBarry Smith   nmax = work[mytid];
13078b31e54SBarry Smith   PETSCFREE(work);
1311eb62cbbSBarry Smith 
1321eb62cbbSBarry Smith   /* post receives:
1331eb62cbbSBarry Smith        1) each message will consist of ordered pairs
1341eb62cbbSBarry Smith      (global index,value) we store the global index as a double
135d6dfbf8fSBarry Smith      to simplify the message passing.
1361eb62cbbSBarry Smith        2) since we don't know how long each individual message is we
1371eb62cbbSBarry Smith      allocate the largest needed buffer for each receive. Potentially
1381eb62cbbSBarry Smith      this is a lot of wasted space.
1391eb62cbbSBarry Smith 
1401eb62cbbSBarry Smith 
1411eb62cbbSBarry Smith        This could be done better.
1421eb62cbbSBarry Smith   */
14378b31e54SBarry Smith   rvalues = (Scalar *) PETSCMALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));
14478b31e54SBarry Smith   CHKPTRQ(rvalues);
14578b31e54SBarry Smith   recv_waits = (MPI_Request *) PETSCMALLOC((nreceives+1)*sizeof(MPI_Request));
14678b31e54SBarry Smith   CHKPTRQ(recv_waits);
1471eb62cbbSBarry Smith   for ( i=0; i<nreceives; i++ ) {
148c60448a5SBarry Smith     MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,
1491eb62cbbSBarry Smith               comm,recv_waits+i);
1501eb62cbbSBarry Smith   }
1511eb62cbbSBarry Smith 
1521eb62cbbSBarry Smith   /* do sends:
1531eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
1541eb62cbbSBarry Smith          the ith processor
1551eb62cbbSBarry Smith   */
15678b31e54SBarry Smith   svalues = (Scalar *) PETSCMALLOC( 3*(aij->stash.n+1)*sizeof(Scalar) );
15778b31e54SBarry Smith   CHKPTRQ(svalues);
15878b31e54SBarry Smith   send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request));
15978b31e54SBarry Smith   CHKPTRQ(send_waits);
16078b31e54SBarry Smith   starts = (int *) PETSCMALLOC( numtids*sizeof(int) ); CHKPTRQ(starts);
1611eb62cbbSBarry Smith   starts[0] = 0;
1621eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1631eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
1641eb62cbbSBarry Smith     svalues[3*starts[owner[i]]]       = (Scalar)  aij->stash.idx[i];
1651eb62cbbSBarry Smith     svalues[3*starts[owner[i]]+1]     = (Scalar)  aij->stash.idy[i];
1661eb62cbbSBarry Smith     svalues[3*(starts[owner[i]]++)+2] =  aij->stash.array[i];
1671eb62cbbSBarry Smith   }
16878b31e54SBarry Smith   PETSCFREE(owner);
1691eb62cbbSBarry Smith   starts[0] = 0;
1701eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1711eb62cbbSBarry Smith   count = 0;
1721eb62cbbSBarry Smith   for ( i=0; i<numtids; i++ ) {
1731eb62cbbSBarry Smith     if (procs[i]) {
174c60448a5SBarry Smith       MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPIU_SCALAR,i,tag,
1751eb62cbbSBarry Smith                 comm,send_waits+count++);
1761eb62cbbSBarry Smith     }
1771eb62cbbSBarry Smith   }
17878b31e54SBarry Smith   PETSCFREE(starts); PETSCFREE(nprocs);
1791eb62cbbSBarry Smith 
1801eb62cbbSBarry Smith   /* Free cache space */
18178b31e54SBarry Smith   ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr);
1821eb62cbbSBarry Smith 
1831eb62cbbSBarry Smith   aij->svalues    = svalues;    aij->rvalues = rvalues;
1841eb62cbbSBarry Smith   aij->nsends     = nsends;     aij->nrecvs = nreceives;
1851eb62cbbSBarry Smith   aij->send_waits = send_waits; aij->recv_waits = recv_waits;
1861eb62cbbSBarry Smith   aij->rmax       = nmax;
1871eb62cbbSBarry Smith 
1881eb62cbbSBarry Smith   return 0;
1891eb62cbbSBarry Smith }
19044a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat);
1911eb62cbbSBarry Smith 
192fc76ce87SLois Curfman McInnes static int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode)
1931eb62cbbSBarry Smith {
1941eb62cbbSBarry Smith   int        ierr;
19544a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1961eb62cbbSBarry Smith   MPI_Status  *send_status,recv_status;
1976abc6512SBarry Smith   int         imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n;
1982493cbb0SBarry Smith   int         row,col,other_disassembled;
1991eb62cbbSBarry Smith   Scalar      *values,val;
2001eb62cbbSBarry Smith   InsertMode  addv = aij->insertmode;
2011eb62cbbSBarry Smith 
2021eb62cbbSBarry Smith   /*  wait on receives */
2031eb62cbbSBarry Smith   while (count) {
204d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);
2051eb62cbbSBarry Smith     /* unpack receives into our local space */
206d6dfbf8fSBarry Smith     values = aij->rvalues + 3*imdex*aij->rmax;
207c60448a5SBarry Smith     MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
2081eb62cbbSBarry Smith     n = n/3;
2091eb62cbbSBarry Smith     for ( i=0; i<n; i++ ) {
2101eb62cbbSBarry Smith       row = (int) PETSCREAL(values[3*i]) - aij->rstart;
2111eb62cbbSBarry Smith       col = (int) PETSCREAL(values[3*i+1]);
2121eb62cbbSBarry Smith       val = values[3*i+2];
2131eb62cbbSBarry Smith       if (col >= aij->cstart && col < aij->cend) {
2141eb62cbbSBarry Smith           col -= aij->cstart;
2151eb62cbbSBarry Smith         MatSetValues(aij->A,1,&row,1,&col,&val,addv);
2161eb62cbbSBarry Smith       }
2171eb62cbbSBarry Smith       else {
218d6dfbf8fSBarry Smith         if (aij->assembled) {
219b7c46309SBarry Smith           if (!aij->colmap) {ierr = CreateColmap_Private(mat);CHKERRQ(ierr);}
2209e25ed09SBarry Smith           col = aij->colmap[col] - 1;
221b7c46309SBarry Smith           if (col < 0  && !((Mat_AIJ*)(aij->A->data))->nonew) {
2222493cbb0SBarry Smith             ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
2232493cbb0SBarry Smith             col = (int) PETSCREAL(values[3*i+1]);
224d6dfbf8fSBarry Smith           }
2259e25ed09SBarry Smith         }
2261eb62cbbSBarry Smith         MatSetValues(aij->B,1,&row,1,&col,&val,addv);
2271eb62cbbSBarry Smith       }
2281eb62cbbSBarry Smith     }
2291eb62cbbSBarry Smith     count--;
2301eb62cbbSBarry Smith   }
23178b31e54SBarry Smith   PETSCFREE(aij->recv_waits); PETSCFREE(aij->rvalues);
2321eb62cbbSBarry Smith 
2331eb62cbbSBarry Smith   /* wait on sends */
2341eb62cbbSBarry Smith   if (aij->nsends) {
23578b31e54SBarry Smith     send_status = (MPI_Status *) PETSCMALLOC( aij->nsends*sizeof(MPI_Status) );
23678b31e54SBarry Smith     CHKPTRQ(send_status);
2371eb62cbbSBarry Smith     MPI_Waitall(aij->nsends,aij->send_waits,send_status);
23878b31e54SBarry Smith     PETSCFREE(send_status);
2391eb62cbbSBarry Smith   }
24078b31e54SBarry Smith   PETSCFREE(aij->send_waits); PETSCFREE(aij->svalues);
2411eb62cbbSBarry Smith 
24207a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
24378b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr);
24478b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr);
2451eb62cbbSBarry Smith 
2462493cbb0SBarry Smith   /* determine if any processor has disassembled, if so we must
2472493cbb0SBarry Smith      also disassemble ourselfs, in order that we may reassemble. */
2482493cbb0SBarry Smith   MPI_Allreduce((void *) &aij->assembled,(void *) &other_disassembled,1,
2492493cbb0SBarry Smith                  MPI_INT,MPI_PROD,mat->comm);
2502493cbb0SBarry Smith   if (aij->assembled && !other_disassembled) {
2512493cbb0SBarry Smith     ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
2522493cbb0SBarry Smith   }
2532493cbb0SBarry Smith 
2545e42470aSBarry Smith   if (!aij->assembled && mode == FINAL_ASSEMBLY) {
25578b31e54SBarry Smith     ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr);
256d6dfbf8fSBarry Smith     aij->assembled = 1;
2575e42470aSBarry Smith   }
25878b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr);
25978b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr);
2605e42470aSBarry Smith 
2618a729477SBarry Smith   return 0;
2628a729477SBarry Smith }
2638a729477SBarry Smith 
2642ee70a88SLois Curfman McInnes static int MatZeroEntries_MPIAIJ(Mat A)
2651eb62cbbSBarry Smith {
26644a69424SLois Curfman McInnes   Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data;
267dbd7a890SLois Curfman McInnes   int ierr;
26878b31e54SBarry Smith   ierr = MatZeroEntries(l->A); CHKERRQ(ierr);
26978b31e54SBarry Smith   ierr = MatZeroEntries(l->B); CHKERRQ(ierr);
2701eb62cbbSBarry Smith   return 0;
2711eb62cbbSBarry Smith }
2721eb62cbbSBarry Smith 
2731eb62cbbSBarry Smith /* again this uses the same basic stratagy as in the assembly and
2741eb62cbbSBarry Smith    scatter create routines, we should try to do it systemamatically
2751eb62cbbSBarry Smith    if we can figure out the proper level of generality. */
2761eb62cbbSBarry Smith 
2771eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the
2781eb62cbbSBarry Smith    matrix is square and the column and row owerships are identical.
2791eb62cbbSBarry Smith    This is a BUG. The only way to fix it seems to be to access
2801eb62cbbSBarry Smith    aij->A and aij->B directly and not through the MatZeroRows()
2811eb62cbbSBarry Smith    routine.
2821eb62cbbSBarry Smith */
28344a69424SLois Curfman McInnes static int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag)
2841eb62cbbSBarry Smith {
28544a69424SLois Curfman McInnes   Mat_MPIAIJ     *l = (Mat_MPIAIJ *) A->data;
2861eb62cbbSBarry Smith   int            i,ierr,N, *rows,*owners = l->rowners,numtids = l->numtids;
2876abc6512SBarry Smith   int            *procs,*nprocs,j,found,idx,nsends,*work;
2881eb62cbbSBarry Smith   int            nmax,*svalues,*starts,*owner,nrecvs,mytid = l->mytid;
2895392566eSBarry Smith   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
290d6dfbf8fSBarry Smith   int            *lens,imdex,*lrows,*values;
291d6dfbf8fSBarry Smith   MPI_Comm       comm = A->comm;
2921eb62cbbSBarry Smith   MPI_Request    *send_waits,*recv_waits;
2931eb62cbbSBarry Smith   MPI_Status     recv_status,*send_status;
2941eb62cbbSBarry Smith   IS             istmp;
2951eb62cbbSBarry Smith 
296c8f67822SLois Curfman McInnes   if (!l->assembled)
29778b31e54SBarry Smith     SETERRQ(1,"MatZeroRows_MPIAIJ:Must assemble matrix first");
29878b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
29978b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
3001eb62cbbSBarry Smith 
3011eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
30278b31e54SBarry Smith   nprocs = (int *) PETSCMALLOC( 2*numtids*sizeof(int) ); CHKPTRQ(nprocs);
30378b31e54SBarry Smith   PETSCMEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids;
30478b31e54SBarry Smith   owner = (int *) PETSCMALLOC((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
3051eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
3061eb62cbbSBarry Smith     idx = rows[i];
3071eb62cbbSBarry Smith     found = 0;
3081eb62cbbSBarry Smith     for ( j=0; j<numtids; j++ ) {
3091eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
3101eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
3111eb62cbbSBarry Smith       }
3121eb62cbbSBarry Smith     }
313bbb6d6a8SBarry Smith     if (!found) SETERRQ(1,"MatZeroRows_MPIAIJ:Index out of range");
3141eb62cbbSBarry Smith   }
3151eb62cbbSBarry Smith   nsends = 0;  for ( i=0; i<numtids; i++ ) { nsends += procs[i];}
3161eb62cbbSBarry Smith 
3171eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
31878b31e54SBarry Smith   work = (int *) PETSCMALLOC( numtids*sizeof(int) ); CHKPTRQ(work);
3191eb62cbbSBarry Smith   MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm);
3201eb62cbbSBarry Smith   nrecvs = work[mytid];
3211eb62cbbSBarry Smith   MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm);
3221eb62cbbSBarry Smith   nmax = work[mytid];
32378b31e54SBarry Smith   PETSCFREE(work);
3241eb62cbbSBarry Smith 
3251eb62cbbSBarry Smith   /* post receives:   */
32678b31e54SBarry Smith   rvalues = (int *) PETSCMALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
32778b31e54SBarry Smith   CHKPTRQ(rvalues);
32878b31e54SBarry Smith   recv_waits = (MPI_Request *) PETSCMALLOC((nrecvs+1)*sizeof(MPI_Request));
32978b31e54SBarry Smith   CHKPTRQ(recv_waits);
3301eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
3311eb62cbbSBarry Smith     MPI_Irecv((void *)(rvalues+nmax*i),nmax,MPI_INT,MPI_ANY_SOURCE,tag,
3321eb62cbbSBarry Smith               comm,recv_waits+i);
3331eb62cbbSBarry Smith   }
3341eb62cbbSBarry Smith 
3351eb62cbbSBarry Smith   /* do sends:
3361eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
3371eb62cbbSBarry Smith          the ith processor
3381eb62cbbSBarry Smith   */
33978b31e54SBarry Smith   svalues = (int *) PETSCMALLOC( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
34078b31e54SBarry Smith   send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request));
34178b31e54SBarry Smith   CHKPTRQ(send_waits);
34278b31e54SBarry Smith   starts = (int *) PETSCMALLOC( (numtids+1)*sizeof(int) ); CHKPTRQ(starts);
3431eb62cbbSBarry Smith   starts[0] = 0;
3441eb62cbbSBarry Smith   for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3451eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
3461eb62cbbSBarry Smith     svalues[starts[owner[i]]++] = rows[i];
3471eb62cbbSBarry Smith   }
3481eb62cbbSBarry Smith   ISRestoreIndices(is,&rows);
3491eb62cbbSBarry Smith 
3501eb62cbbSBarry Smith   starts[0] = 0;
3511eb62cbbSBarry Smith   for ( i=1; i<numtids+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3521eb62cbbSBarry Smith   count = 0;
3531eb62cbbSBarry Smith   for ( i=0; i<numtids; i++ ) {
3541eb62cbbSBarry Smith     if (procs[i]) {
3551eb62cbbSBarry Smith       MPI_Isend((void*)(svalues+starts[i]),nprocs[i],MPI_INT,i,tag,
3561eb62cbbSBarry Smith                 comm,send_waits+count++);
3571eb62cbbSBarry Smith     }
3581eb62cbbSBarry Smith   }
35978b31e54SBarry Smith   PETSCFREE(starts);
3601eb62cbbSBarry Smith 
3611eb62cbbSBarry Smith   base = owners[mytid];
3621eb62cbbSBarry Smith 
3631eb62cbbSBarry Smith   /*  wait on receives */
36478b31e54SBarry Smith   lens = (int *) PETSCMALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
3651eb62cbbSBarry Smith   source = lens + nrecvs;
3661eb62cbbSBarry Smith   count = nrecvs; slen = 0;
3671eb62cbbSBarry Smith   while (count) {
368d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
3691eb62cbbSBarry Smith     /* unpack receives into our local space */
3701eb62cbbSBarry Smith     MPI_Get_count(&recv_status,MPI_INT,&n);
371d6dfbf8fSBarry Smith     source[imdex]  = recv_status.MPI_SOURCE;
372d6dfbf8fSBarry Smith     lens[imdex]  = n;
3731eb62cbbSBarry Smith     slen += n;
3741eb62cbbSBarry Smith     count--;
3751eb62cbbSBarry Smith   }
37678b31e54SBarry Smith   PETSCFREE(recv_waits);
3771eb62cbbSBarry Smith 
3781eb62cbbSBarry Smith   /* move the data into the send scatter */
37978b31e54SBarry Smith   lrows = (int *) PETSCMALLOC( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
3801eb62cbbSBarry Smith   count = 0;
3811eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
3821eb62cbbSBarry Smith     values = rvalues + i*nmax;
3831eb62cbbSBarry Smith     for ( j=0; j<lens[i]; j++ ) {
3841eb62cbbSBarry Smith       lrows[count++] = values[j] - base;
3851eb62cbbSBarry Smith     }
3861eb62cbbSBarry Smith   }
38778b31e54SBarry Smith   PETSCFREE(rvalues); PETSCFREE(lens);
38878b31e54SBarry Smith   PETSCFREE(owner); PETSCFREE(nprocs);
3891eb62cbbSBarry Smith 
3901eb62cbbSBarry Smith   /* actually zap the local rows */
391*fafbff53SBarry Smith   ierr = ISCreateSeq(MPI_COMM_SELF,slen,lrows,&istmp);
392464493b3SBarry Smith   PLogObjectParent(A,istmp);
39378b31e54SBarry Smith   CHKERRQ(ierr);  PETSCFREE(lrows);
39478b31e54SBarry Smith   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
39578b31e54SBarry Smith   ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr);
39678b31e54SBarry Smith   ierr = ISDestroy(istmp); CHKERRQ(ierr);
3971eb62cbbSBarry Smith 
3981eb62cbbSBarry Smith   /* wait on sends */
3991eb62cbbSBarry Smith   if (nsends) {
40078b31e54SBarry Smith     send_status = (MPI_Status *) PETSCMALLOC( nsends*sizeof(MPI_Status) );
40178b31e54SBarry Smith     CHKPTRQ(send_status);
4021eb62cbbSBarry Smith     MPI_Waitall(nsends,send_waits,send_status);
40378b31e54SBarry Smith     PETSCFREE(send_status);
4041eb62cbbSBarry Smith   }
40578b31e54SBarry Smith   PETSCFREE(send_waits); PETSCFREE(svalues);
4061eb62cbbSBarry Smith 
4071eb62cbbSBarry Smith 
4081eb62cbbSBarry Smith   return 0;
4091eb62cbbSBarry Smith }
4101eb62cbbSBarry Smith 
41144a69424SLois Curfman McInnes static int MatMult_MPIAIJ(Mat aijin,Vec xx,Vec yy)
4121eb62cbbSBarry Smith {
41344a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
4141eb62cbbSBarry Smith   int        ierr;
41578b31e54SBarry Smith   if (!aij->assembled) SETERRQ(1,"MatMult_MPIAIJ:must assemble matrix first");
41606be10caSBarry Smith   ierr = VecScatterBegin(xx,aij->lvec,INSERTVALUES,SCATTERALL,aij->Mvctx);
41778b31e54SBarry Smith   CHKERRQ(ierr);
41878b31e54SBarry Smith   ierr = MatMult(aij->A,xx,yy); CHKERRQ(ierr);
41906be10caSBarry Smith   ierr = VecScatterEnd(xx,aij->lvec,INSERTVALUES,SCATTERALL,aij->Mvctx);
42078b31e54SBarry Smith   CHKERRQ(ierr);
42178b31e54SBarry Smith   ierr = MatMultAdd(aij->B,aij->lvec,yy,yy); CHKERRQ(ierr);
4221eb62cbbSBarry Smith   return 0;
4231eb62cbbSBarry Smith }
4241eb62cbbSBarry Smith 
42544a69424SLois Curfman McInnes static int MatMultAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz)
426da3a660dSBarry Smith {
42744a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
428da3a660dSBarry Smith   int        ierr;
42978b31e54SBarry Smith   if (!aij->assembled) SETERRQ(1,"MatMult_MPIAIJ:must assemble matrix first");
43006be10caSBarry Smith   ierr = VecScatterBegin(xx,aij->lvec,INSERTVALUES,SCATTERALL,aij->Mvctx);
43178b31e54SBarry Smith   CHKERRQ(ierr);
43278b31e54SBarry Smith   ierr = MatMultAdd(aij->A,xx,yy,zz); CHKERRQ(ierr);
43306be10caSBarry Smith   ierr = VecScatterEnd(xx,aij->lvec,INSERTVALUES,SCATTERALL,aij->Mvctx);
43478b31e54SBarry Smith   CHKERRQ(ierr);
43578b31e54SBarry Smith   ierr = MatMultAdd(aij->B,aij->lvec,zz,zz); CHKERRQ(ierr);
436da3a660dSBarry Smith   return 0;
437da3a660dSBarry Smith }
438da3a660dSBarry Smith 
43944a69424SLois Curfman McInnes static int MatMultTrans_MPIAIJ(Mat aijin,Vec xx,Vec yy)
440da3a660dSBarry Smith {
44144a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
442da3a660dSBarry Smith   int        ierr;
443da3a660dSBarry Smith 
44444a69424SLois Curfman McInnes   if (!aij->assembled)
44578b31e54SBarry Smith     SETERRQ(1,"MatMulTrans_MPIAIJ:must assemble matrix first");
446da3a660dSBarry Smith   /* do nondiagonal part */
44778b31e54SBarry Smith   ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERRQ(ierr);
448da3a660dSBarry Smith   /* send it on its way */
44906be10caSBarry Smith   ierr = VecScatterBegin(aij->lvec,yy,ADDVALUES,
45078b31e54SBarry Smith            (ScatterMode)(SCATTERALL|SCATTERREVERSE),aij->Mvctx); CHKERRQ(ierr);
451da3a660dSBarry Smith   /* do local part */
45278b31e54SBarry Smith   ierr = MatMultTrans(aij->A,xx,yy); 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. */
45606be10caSBarry Smith   ierr = VecScatterEnd(aij->lvec,yy,ADDVALUES,
45778b31e54SBarry Smith          (ScatterMode)(SCATTERALL|SCATTERREVERSE),aij->Mvctx); CHKERRQ(ierr);
458da3a660dSBarry Smith   return 0;
459da3a660dSBarry Smith }
460da3a660dSBarry Smith 
46144a69424SLois Curfman McInnes static int MatMultTransAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz)
462da3a660dSBarry Smith {
46344a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
464da3a660dSBarry Smith   int        ierr;
465da3a660dSBarry Smith 
46644a69424SLois Curfman McInnes   if (!aij->assembled)
46778b31e54SBarry Smith     SETERRQ(1,"MatMulTransAdd_MPIAIJ:must assemble matrix first");
468da3a660dSBarry Smith   /* do nondiagonal part */
46978b31e54SBarry Smith   ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERRQ(ierr);
470da3a660dSBarry Smith   /* send it on its way */
47106be10caSBarry Smith   ierr = VecScatterBegin(aij->lvec,zz,ADDVALUES,
47278b31e54SBarry Smith          (ScatterMode)(SCATTERALL|SCATTERREVERSE),aij->Mvctx); CHKERRQ(ierr);
473da3a660dSBarry Smith   /* do local part */
47478b31e54SBarry Smith   ierr = MatMultTransAdd(aij->A,xx,yy,zz); CHKERRQ(ierr);
475da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
476da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
477da3a660dSBarry Smith   /* but is not perhaps always true. */
47806be10caSBarry Smith   ierr = VecScatterEnd(aij->lvec,zz,ADDVALUES,
47978b31e54SBarry Smith        (ScatterMode)(SCATTERALL|SCATTERREVERSE),aij->Mvctx); CHKERRQ(ierr);
480da3a660dSBarry Smith   return 0;
481da3a660dSBarry Smith }
482da3a660dSBarry Smith 
4831eb62cbbSBarry Smith /*
4841eb62cbbSBarry Smith   This only works correctly for square matrices where the subblock A->A is the
4851eb62cbbSBarry Smith    diagonal block
4861eb62cbbSBarry Smith */
487d1710a03SLois Curfman McInnes static int MatGetDiagonal_MPIAIJ(Mat Ain,Vec v)
4881eb62cbbSBarry Smith {
48944a69424SLois Curfman McInnes   Mat_MPIAIJ *A = (Mat_MPIAIJ *) Ain->data;
49078b31e54SBarry Smith   if (!A->assembled) SETERRQ(1,"MatGetDiag_MPIAIJ:must assemble matrix first");
4911eb62cbbSBarry Smith   return MatGetDiagonal(A->A,v);
4921eb62cbbSBarry Smith }
4931eb62cbbSBarry Smith 
49444a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj)
4951eb62cbbSBarry Smith {
4961eb62cbbSBarry Smith   Mat        mat = (Mat) obj;
49744a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
4981eb62cbbSBarry Smith   int        ierr;
499a5a9c739SBarry Smith #if defined(PETSC_LOG)
5006e35fa57SLois Curfman McInnes   PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N);
501a5a9c739SBarry Smith #endif
50278b31e54SBarry Smith   PETSCFREE(aij->rowners);
50378b31e54SBarry Smith   ierr = MatDestroy(aij->A); CHKERRQ(ierr);
50478b31e54SBarry Smith   ierr = MatDestroy(aij->B); CHKERRQ(ierr);
50578b31e54SBarry Smith   if (aij->colmap) PETSCFREE(aij->colmap);
50678b31e54SBarry Smith   if (aij->garray) PETSCFREE(aij->garray);
5071eb62cbbSBarry Smith   if (aij->lvec) VecDestroy(aij->lvec);
5081eb62cbbSBarry Smith   if (aij->Mvctx) VecScatterCtxDestroy(aij->Mvctx);
50978b31e54SBarry Smith   PETSCFREE(aij);
510a5a9c739SBarry Smith   PLogObjectDestroy(mat);
511a5a9c739SBarry Smith   PETSCHEADERDESTROY(mat);
5121eb62cbbSBarry Smith   return 0;
5131eb62cbbSBarry Smith }
5147857610eSBarry Smith #include "draw.h"
515b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
516ee50ffe9SBarry Smith 
51744a69424SLois Curfman McInnes static int MatView_MPIAIJ(PetscObject obj,Viewer viewer)
5181eb62cbbSBarry Smith {
5191eb62cbbSBarry Smith   Mat         mat = (Mat) obj;
52044a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
521f72585f2SLois Curfman McInnes   int         ierr, format;
522d13ab20cSBarry Smith   PetscObject vobj = (PetscObject) viewer;
5231eb62cbbSBarry Smith 
52478b31e54SBarry Smith   if (!aij->assembled) SETERRQ(1,"MatView_MPIAIJ:must assemble matrix first");
525154123eaSLois Curfman McInnes   if (!viewer) { /* so that viewers may be used from debuggers */
52621b0d8fbSLois Curfman McInnes     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
527154123eaSLois Curfman McInnes   }
528ab254492SBarry Smith   if (vobj->cookie == DRAW_COOKIE && vobj->type == NULLWINDOW) return 0;
529d636dbe3SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
530f72585f2SLois Curfman McInnes   if (vobj->cookie == VIEWER_COOKIE && format == FILE_FORMAT_INFO &&
5319b94acceSBarry Smith      (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER)) {
532f72585f2SLois Curfman McInnes    /* do nothing for now */
533f72585f2SLois Curfman McInnes   }
5349b94acceSBarry Smith   else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILE_VIEWER) {
535d636dbe3SBarry Smith     FILE *fd;
536d636dbe3SBarry Smith     ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
5377f813858SBarry Smith     MPIU_Seq_begin(mat->comm,1);
538d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
5391eb62cbbSBarry Smith            aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
5401eb62cbbSBarry Smith            aij->cend);
54178b31e54SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
54278b31e54SBarry Smith     ierr = MatView(aij->B,viewer); CHKERRQ(ierr);
543d13ab20cSBarry Smith     fflush(fd);
5447f813858SBarry Smith     MPIU_Seq_end(mat->comm,1);
545d13ab20cSBarry Smith   }
5469b94acceSBarry Smith   else if ((vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILES_VIEWER) ||
5477857610eSBarry Smith             vobj->cookie == DRAW_COOKIE) {
54895373324SBarry Smith     int numtids = aij->numtids, mytid = aij->mytid;
54995373324SBarry Smith     if (numtids == 1) {
55078b31e54SBarry Smith       ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
55195373324SBarry Smith     }
55295373324SBarry Smith     else {
55395373324SBarry Smith       /* assemble the entire matrix onto first processor. */
55495373324SBarry Smith       Mat     A;
5552ee70a88SLois Curfman McInnes       Mat_AIJ *Aloc;
5562eb8c8abSBarry Smith       int     M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
55795373324SBarry Smith       Scalar  *a;
5582ee70a88SLois Curfman McInnes 
55995373324SBarry Smith       if (!mytid) {
56095373324SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,0,0,0,&A);
56195373324SBarry Smith       }
56295373324SBarry Smith       else {
56395373324SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,0,0,0,&A);
56495373324SBarry Smith       }
565464493b3SBarry Smith       PLogObjectParent(mat,A);
56678b31e54SBarry Smith       CHKERRQ(ierr);
56795373324SBarry Smith 
56895373324SBarry Smith       /* copy over the A part */
5692ee70a88SLois Curfman McInnes       Aloc = (Mat_AIJ*) aij->A->data;
5702ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
57195373324SBarry Smith       row = aij->rstart;
572b7c46309SBarry Smith       for ( i=0; i<ai[m]-1; i++ ) {aj[i] += aij->cstart - 1;}
57395373324SBarry Smith       for ( i=0; i<m; i++ ) {
57407a0e7adSLois Curfman McInnes         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERTVALUES);
57578b31e54SBarry Smith         CHKERRQ(ierr);
57695373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
57795373324SBarry Smith       }
5782ee70a88SLois Curfman McInnes       aj = Aloc->j;
579b7c46309SBarry Smith       for ( i=0; i<ai[m]-1; i++ ) {aj[i] -= aij->cstart - 1;}
58095373324SBarry Smith 
58195373324SBarry Smith       /* copy over the B part */
5822ee70a88SLois Curfman McInnes       Aloc = (Mat_AIJ*) aij->B->data;
5832ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
58495373324SBarry Smith       row = aij->rstart;
58578b31e54SBarry Smith       ct = cols = (int *) PETSCMALLOC( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
586b7c46309SBarry Smith       for ( i=0; i<ai[m]-1; i++ ) {cols[i] = aij->garray[aj[i]-1];}
58795373324SBarry Smith       for ( i=0; i<m; i++ ) {
58807a0e7adSLois Curfman McInnes         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERTVALUES);
58978b31e54SBarry Smith         CHKERRQ(ierr);
59095373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
59195373324SBarry Smith       }
59278b31e54SBarry Smith       PETSCFREE(ct);
59378b31e54SBarry Smith       ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
59478b31e54SBarry Smith       ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
59595373324SBarry Smith       if (!mytid) {
59678b31e54SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr);
59795373324SBarry Smith       }
59878b31e54SBarry Smith       ierr = MatDestroy(A); CHKERRQ(ierr);
59995373324SBarry Smith     }
60095373324SBarry Smith   }
6011eb62cbbSBarry Smith   return 0;
6021eb62cbbSBarry Smith }
6031eb62cbbSBarry Smith 
60480b4ade8SLois Curfman McInnes extern int MatMarkDiag_AIJ(Mat);
6051eb62cbbSBarry Smith /*
6061eb62cbbSBarry Smith     This has to provide several versions.
6071eb62cbbSBarry Smith 
6081eb62cbbSBarry Smith      1) per sequential
6091eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
6101eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
611d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
6121eb62cbbSBarry Smith */
613fc76ce87SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
614fc76ce87SLois Curfman McInnes                            double shift,int its,Vec xx)
6158a729477SBarry Smith {
61644a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
617d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
61844a69424SLois Curfman McInnes   Mat_AIJ    *A = (Mat_AIJ *) AA->data, *B = (Mat_AIJ *)BB->data;
6196abc6512SBarry Smith   Scalar     zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts;
6206abc6512SBarry Smith   int        ierr,*idx, *diag;
6216abc6512SBarry Smith   int        n = mat->n, m = mat->m, i;
622da3a660dSBarry Smith   Vec        tt;
6238a729477SBarry Smith 
62478b31e54SBarry Smith   if (!mat->assembled) SETERRQ(1,"MatRelax_MPIAIJ:must assemble matrix first");
6251eb62cbbSBarry Smith 
626d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
627d6dfbf8fSBarry Smith   xs = x -1; /* shift by one for index start of 1 */
628da3a660dSBarry Smith   ls--;
62980b4ade8SLois Curfman McInnes   if (!A->diag) {if ((ierr = MatMarkDiag_AIJ(AA))) return ierr;}
630d6dfbf8fSBarry Smith   diag = A->diag;
631acb40c82SBarry Smith   if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) {
63242f8db3fSLois Curfman McInnes     SETERRQ(1,"MatRelax_MPIAIJ:Option not yet supported");
633acb40c82SBarry Smith   }
634da3a660dSBarry Smith   if (flag & SOR_EISENSTAT) {
635da3a660dSBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
636da3a660dSBarry Smith     U is upper triangular, E is diagonal; This routine applies
637da3a660dSBarry Smith 
638da3a660dSBarry Smith             (L + E)^{-1} A (U + E)^{-1}
639da3a660dSBarry Smith 
640da3a660dSBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
641da3a660dSBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
642da3a660dSBarry Smith     is the relaxation factor.
643da3a660dSBarry Smith     */
64478b31e54SBarry Smith     ierr = VecDuplicate(xx,&tt); CHKERRQ(ierr);
645464493b3SBarry Smith     PLogObjectParent(matin,tt);
646da3a660dSBarry Smith     VecGetArray(tt,&t);
647da3a660dSBarry Smith     scale = (2.0/omega) - 1.0;
648da3a660dSBarry Smith     /*  x = (E + U)^{-1} b */
649da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
65006be10caSBarry Smith     ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
65178b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
652da3a660dSBarry Smith     for ( i=m-1; i>-1; i-- ) {
653da3a660dSBarry Smith       n    = A->i[i+1] - diag[i] - 1;
654da3a660dSBarry Smith       idx  = A->j + diag[i];
655da3a660dSBarry Smith       v    = A->a + diag[i];
656da3a660dSBarry Smith       sum  = b[i];
657da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
658da3a660dSBarry Smith       d    = shift + A->a[diag[i]-1];
659da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
660da3a660dSBarry Smith       idx  = B->j + B->i[i] - 1;
661da3a660dSBarry Smith       v    = B->a + B->i[i] - 1;
662da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
663da3a660dSBarry Smith       x[i] = omega*(sum/d);
664da3a660dSBarry Smith     }
665edd2f0e1SBarry Smith     ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
66678b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
667da3a660dSBarry Smith 
668da3a660dSBarry Smith     /*  t = b - (2*E - D)x */
669da3a660dSBarry Smith     v = A->a;
670da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ - 1])*x[i]; }
671da3a660dSBarry Smith 
672da3a660dSBarry Smith     /*  t = (E + L)^{-1}t */
673da3a660dSBarry Smith     ts = t - 1; /* shifted by one for index start of a or mat->j*/
674da3a660dSBarry Smith     diag = A->diag;
675da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
67606be10caSBarry Smith     ierr = VecPipelineBegin(tt,mat->lvec,INSERTVALUES,PIPELINEDOWN,
67778b31e54SBarry Smith                                                  mat->Mvctx); CHKERRQ(ierr);
678da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
679da3a660dSBarry Smith       n    = diag[i] - A->i[i];
680da3a660dSBarry Smith       idx  = A->j + A->i[i] - 1;
681da3a660dSBarry Smith       v    = A->a + A->i[i] - 1;
682da3a660dSBarry Smith       sum  = t[i];
683da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
684da3a660dSBarry Smith       d    = shift + A->a[diag[i]-1];
685da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
686da3a660dSBarry Smith       idx  = B->j + B->i[i] - 1;
687da3a660dSBarry Smith       v    = B->a + B->i[i] - 1;
688da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
689da3a660dSBarry Smith       t[i] = omega*(sum/d);
690da3a660dSBarry Smith     }
691edd2f0e1SBarry Smith     ierr = VecPipelineEnd(tt,mat->lvec,INSERTVALUES,PIPELINEDOWN,
69278b31e54SBarry Smith                                                     mat->Mvctx); CHKERRQ(ierr);
693da3a660dSBarry Smith     /*  x = x + t */
694da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
695da3a660dSBarry Smith     VecDestroy(tt);
696da3a660dSBarry Smith     return 0;
697da3a660dSBarry Smith   }
698da3a660dSBarry Smith 
6991eb62cbbSBarry Smith 
700d6dfbf8fSBarry Smith   if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){
701da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
702da3a660dSBarry Smith       VecSet(&zero,mat->lvec); VecSet(&zero,xx);
703da3a660dSBarry Smith     }
704da3a660dSBarry Smith     else {
70506be10caSBarry Smith       ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
70678b31e54SBarry Smith       CHKERRQ(ierr);
70706be10caSBarry Smith       ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
70878b31e54SBarry Smith       CHKERRQ(ierr);
709da3a660dSBarry Smith     }
710d6dfbf8fSBarry Smith     while (its--) {
711d6dfbf8fSBarry Smith       /* go down through the rows */
71206be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
71378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
714d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
715d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
716d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
717d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
718d6dfbf8fSBarry Smith         sum  = b[i];
719d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
720d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
721d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
722d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
723d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
724d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
725d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
726d6dfbf8fSBarry Smith       }
727edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
72878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
729d6dfbf8fSBarry Smith       /* come up through the rows */
73006be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
73178b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
732d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; 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       }
745edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
74678b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
747d6dfbf8fSBarry Smith     }
748d6dfbf8fSBarry Smith   }
749d6dfbf8fSBarry Smith   else if (flag & SOR_FORWARD_SWEEP){
750da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
751da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
75206be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
75378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
754da3a660dSBarry Smith       for ( i=0; i<m; i++ ) {
755da3a660dSBarry Smith         n    = diag[i] - A->i[i];
756da3a660dSBarry Smith         idx  = A->j + A->i[i] - 1;
757da3a660dSBarry Smith         v    = A->a + A->i[i] - 1;
758da3a660dSBarry Smith         sum  = b[i];
759da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
760da3a660dSBarry Smith         d    = shift + A->a[diag[i]-1];
761da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
762da3a660dSBarry Smith         idx  = B->j + B->i[i] - 1;
763da3a660dSBarry Smith         v    = B->a + B->i[i] - 1;
764da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
765da3a660dSBarry Smith         x[i] = omega*(sum/d);
766da3a660dSBarry Smith       }
767edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
76878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
769da3a660dSBarry Smith       its--;
770da3a660dSBarry Smith     }
771d6dfbf8fSBarry Smith     while (its--) {
77206be10caSBarry Smith       ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
77378b31e54SBarry Smith       CHKERRQ(ierr);
77406be10caSBarry Smith       ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
77578b31e54SBarry Smith       CHKERRQ(ierr);
77606be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
77778b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
778d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
779d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
780d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
781d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
782d6dfbf8fSBarry Smith         sum  = b[i];
783d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
784d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
785d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
786d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
787d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
788d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
789d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
790d6dfbf8fSBarry Smith       }
791edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
79278b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
793d6dfbf8fSBarry Smith     }
794d6dfbf8fSBarry Smith   }
795d6dfbf8fSBarry Smith   else if (flag & SOR_BACKWARD_SWEEP){
796da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
797da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
79806be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
79978b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
800da3a660dSBarry Smith       for ( i=m-1; i>-1; i-- ) {
801da3a660dSBarry Smith         n    = A->i[i+1] - diag[i] - 1;
802da3a660dSBarry Smith         idx  = A->j + diag[i];
803da3a660dSBarry Smith         v    = A->a + diag[i];
804da3a660dSBarry Smith         sum  = b[i];
805da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
806da3a660dSBarry Smith         d    = shift + A->a[diag[i]-1];
807da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
808da3a660dSBarry Smith         idx  = B->j + B->i[i] - 1;
809da3a660dSBarry Smith         v    = B->a + B->i[i] - 1;
810da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
811da3a660dSBarry Smith         x[i] = omega*(sum/d);
812da3a660dSBarry Smith       }
813edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
81478b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
815da3a660dSBarry Smith       its--;
816da3a660dSBarry Smith     }
817d6dfbf8fSBarry Smith     while (its--) {
81806be10caSBarry Smith       ierr = VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERDOWN,
81978b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
82006be10caSBarry Smith       ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERDOWN,
82178b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
82206be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
82378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
824d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
825d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
826d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
827d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
828d6dfbf8fSBarry Smith         sum  = b[i];
829d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
830d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
831d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
832d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
833d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
834d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
835d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
836d6dfbf8fSBarry Smith       }
837edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
83878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
839d6dfbf8fSBarry Smith     }
840d6dfbf8fSBarry Smith   }
841d6dfbf8fSBarry Smith   else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
842da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
843da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
844da3a660dSBarry Smith     }
84506be10caSBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
84678b31e54SBarry Smith     CHKERRQ(ierr);
84706be10caSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
84878b31e54SBarry Smith     CHKERRQ(ierr);
849d6dfbf8fSBarry Smith     while (its--) {
850d6dfbf8fSBarry Smith       /* go down through the rows */
851d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
852d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
853d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
854d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
855d6dfbf8fSBarry Smith         sum  = b[i];
856d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
857d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
858d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
859d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
860d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
861d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
862d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
863d6dfbf8fSBarry Smith       }
864d6dfbf8fSBarry Smith       /* come up through the rows */
865d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
866d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
867d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
868d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
869d6dfbf8fSBarry Smith         sum  = b[i];
870d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
871d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
872d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
873d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
874d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
875d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
876d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
877d6dfbf8fSBarry Smith       }
878d6dfbf8fSBarry Smith     }
879d6dfbf8fSBarry Smith   }
880d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_FORWARD_SWEEP){
881da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
882da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
883da3a660dSBarry Smith     }
88406be10caSBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
88578b31e54SBarry Smith     CHKERRQ(ierr);
88606be10caSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
88778b31e54SBarry Smith     CHKERRQ(ierr);
888d6dfbf8fSBarry Smith     while (its--) {
889d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
890d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
891d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
892d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
893d6dfbf8fSBarry Smith         sum  = b[i];
894d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
895d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
896d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
897d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
898d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
899d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
900d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
901d6dfbf8fSBarry Smith       }
902d6dfbf8fSBarry Smith     }
903d6dfbf8fSBarry Smith   }
904d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_BACKWARD_SWEEP){
905da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
906da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
907da3a660dSBarry Smith     }
90806be10caSBarry Smith     ierr = VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERALL,
90978b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
91006be10caSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERALL,
91178b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
912d6dfbf8fSBarry Smith     while (its--) {
913d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
914d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
915d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
916d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
917d6dfbf8fSBarry Smith         sum  = b[i];
918d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
919d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
920d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
921d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
922d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
923d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
924d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
925d6dfbf8fSBarry Smith       }
926d6dfbf8fSBarry Smith     }
927d6dfbf8fSBarry Smith   }
9288a729477SBarry Smith   return 0;
9298a729477SBarry Smith }
930a66be287SLois Curfman McInnes 
931fc76ce87SLois Curfman McInnes static int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,int *nz,
932fc76ce87SLois Curfman McInnes                              int *nzalloc,int *mem)
933a66be287SLois Curfman McInnes {
934a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
935a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
936a66be287SLois Curfman McInnes   int        ierr, isend[3], irecv[3], nzA, nzallocA, memA;
937a66be287SLois Curfman McInnes 
93878b31e54SBarry Smith   ierr = MatGetInfo(A,MAT_LOCAL,&nzA,&nzallocA,&memA); CHKERRQ(ierr);
93978b31e54SBarry Smith   ierr = MatGetInfo(B,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr);
940a66be287SLois Curfman McInnes   isend[0] += nzA; isend[1] += nzallocA; isend[2] += memA;
941a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
942a66be287SLois Curfman McInnes     *nz = isend[0]; *nzalloc = isend[1]; *mem = isend[2];
943a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
944a66be287SLois Curfman McInnes     MPI_Allreduce((void *) isend,(void *) irecv,3,MPI_INT,MPI_MAX,matin->comm);
945a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
946a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
947a66be287SLois Curfman McInnes     MPI_Allreduce((void *) isend,(void *) irecv,3,MPI_INT,MPI_SUM,matin->comm);
948a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
949a66be287SLois Curfman McInnes   }
950a66be287SLois Curfman McInnes   return 0;
951a66be287SLois Curfman McInnes }
952a66be287SLois Curfman McInnes 
953299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*);
954299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*);
955299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double);
956299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *);
957299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec);
958299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec);
959299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec);
960299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec);
961299609e3SLois Curfman McInnes 
962fc76ce87SLois Curfman McInnes static int MatSetOption_MPIAIJ(Mat aijin,MatOption op)
963c74985f6SBarry Smith {
96444a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
965c74985f6SBarry Smith 
966c74985f6SBarry Smith   if      (op == NO_NEW_NONZERO_LOCATIONS)  {
967c74985f6SBarry Smith     MatSetOption(aij->A,op);
968c74985f6SBarry Smith     MatSetOption(aij->B,op);
969c74985f6SBarry Smith   }
970c74985f6SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) {
971c74985f6SBarry Smith     MatSetOption(aij->A,op);
972c74985f6SBarry Smith     MatSetOption(aij->B,op);
973c74985f6SBarry Smith   }
974bbb6d6a8SBarry Smith   else if (op == COLUMN_ORIENTED)
975bbb6d6a8SBarry Smith     SETERRQ(1,"MatSetOption_MPIAIJ:Column oriented not supported");
976c74985f6SBarry Smith   return 0;
977c74985f6SBarry Smith }
978c74985f6SBarry Smith 
979d1710a03SLois Curfman McInnes static int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
980c74985f6SBarry Smith {
98144a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
982c74985f6SBarry Smith   *m = mat->M; *n = mat->N;
983c74985f6SBarry Smith   return 0;
984c74985f6SBarry Smith }
985c74985f6SBarry Smith 
986d1710a03SLois Curfman McInnes static int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
987c74985f6SBarry Smith {
98844a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
989b7c46309SBarry Smith   *m = mat->m; *n = mat->N;
990c74985f6SBarry Smith   return 0;
991c74985f6SBarry Smith }
992c74985f6SBarry Smith 
993d1710a03SLois Curfman McInnes static int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
994c74985f6SBarry Smith {
99544a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
996c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
997c74985f6SBarry Smith   return 0;
998c74985f6SBarry Smith }
999c74985f6SBarry Smith 
100039e00950SLois Curfman McInnes static int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
100139e00950SLois Curfman McInnes {
1002154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1003154123eaSLois Curfman McInnes   Scalar     *vworkA, *vworkB, **pvA, **pvB;
1004154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
1005154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
100639e00950SLois Curfman McInnes 
100739e00950SLois Curfman McInnes   if (row < rstart || row >= rend)
1008bbb6d6a8SBarry Smith     SETERRQ(1,"MatGetRow_MPIAIJ:Currently you can get only local rows")
1009abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
101039e00950SLois Curfman McInnes 
1011154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1012154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1013154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
101478b31e54SBarry Smith   ierr = MatGetRow(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
101578b31e54SBarry Smith   ierr = MatGetRow(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
1016154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1017154123eaSLois Curfman McInnes 
1018154123eaSLois Curfman McInnes   if (v  || idx) {
1019154123eaSLois Curfman McInnes     if (nztot) {
1020154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
1021299609e3SLois Curfman McInnes       int imark;
102248b35521SBarry Smith       if (mat->assembled) {
1023154123eaSLois Curfman McInnes         for (i=0; i<nzB; i++) cworkB[i] = mat->garray[cworkB[i]];
102448b35521SBarry Smith       }
1025154123eaSLois Curfman McInnes       if (v) {
102678b31e54SBarry Smith         *v = (Scalar *) PETSCMALLOC( (nztot)*sizeof(Scalar) ); CHKPTRQ(*v);
102739e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1028154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*v)[i] = vworkB[i];
1029154123eaSLois Curfman McInnes           else break;
1030154123eaSLois Curfman McInnes         }
1031154123eaSLois Curfman McInnes         imark = i;
1032154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*v)[imark+i] = vworkA[i];
1033299609e3SLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*v)[nzA+i] = vworkB[i];
1034154123eaSLois Curfman McInnes       }
1035154123eaSLois Curfman McInnes       if (idx) {
103678b31e54SBarry Smith         *idx = (int *) PETSCMALLOC( (nztot)*sizeof(int) ); CHKPTRQ(*idx);
1037154123eaSLois Curfman McInnes         for (i=0; i<nzA; i++) cworkA[i] += cstart;
1038154123eaSLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1039154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*idx)[i] = cworkB[i];
1040154123eaSLois Curfman McInnes           else break;
1041154123eaSLois Curfman McInnes         }
1042154123eaSLois Curfman McInnes         imark = i;
1043154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*idx)[imark+i] = cworkA[i];
1044299609e3SLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*idx)[nzA+i] = cworkB[i];
104539e00950SLois Curfman McInnes       }
104639e00950SLois Curfman McInnes     }
104739e00950SLois Curfman McInnes     else {*idx = 0; *v=0;}
1048154123eaSLois Curfman McInnes   }
104939e00950SLois Curfman McInnes   *nz = nztot;
105078b31e54SBarry Smith   ierr = MatRestoreRow(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
105178b31e54SBarry Smith   ierr = MatRestoreRow(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
105239e00950SLois Curfman McInnes   return 0;
105339e00950SLois Curfman McInnes }
105439e00950SLois Curfman McInnes 
105539e00950SLois Curfman McInnes static int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
105639e00950SLois Curfman McInnes {
105778b31e54SBarry Smith   if (idx) PETSCFREE(*idx);
105878b31e54SBarry Smith   if (v) PETSCFREE(*v);
105939e00950SLois Curfman McInnes   return 0;
106039e00950SLois Curfman McInnes }
106139e00950SLois Curfman McInnes 
1062855ac2c5SLois Curfman McInnes static int MatNorm_MPIAIJ(Mat mat,MatNormType type,double *norm)
1063855ac2c5SLois Curfman McInnes {
1064855ac2c5SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1065aac3204fSLois Curfman McInnes   int        ierr, i, j, cstart = aij->cstart;
106604ca555eSLois Curfman McInnes   double     sum = 0.0;
106704ca555eSLois Curfman McInnes   Mat_AIJ    *amat = (Mat_AIJ*) aij->A->data, *bmat = (Mat_AIJ*) aij->B->data;
106804ca555eSLois Curfman McInnes   Scalar     *v;
106904ca555eSLois Curfman McInnes 
107004ca555eSLois Curfman McInnes   if (!aij->assembled)
107104ca555eSLois Curfman McInnes     SETERRQ(1,"MatNorm_MPIAIJ: Cannot compute norm of unassembled matrix");
1072855ac2c5SLois Curfman McInnes   if (aij->numtids == 1) {
107314183eadSLois Curfman McInnes     ierr =  MatNorm(aij->A,type,norm); CHKERRQ(ierr);
107437fa93a5SLois Curfman McInnes   } else {
107504ca555eSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
107604ca555eSLois Curfman McInnes       v = amat->a;
107704ca555eSLois Curfman McInnes       for (i=0; i<amat->nz; i++ ) {
107804ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
107904ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
108004ca555eSLois Curfman McInnes #else
108104ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
108204ca555eSLois Curfman McInnes #endif
108304ca555eSLois Curfman McInnes       }
108404ca555eSLois Curfman McInnes       v = bmat->a;
108504ca555eSLois Curfman McInnes       for (i=0; i<bmat->nz; i++ ) {
108604ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
108704ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
108804ca555eSLois Curfman McInnes #else
108904ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
109004ca555eSLois Curfman McInnes #endif
109104ca555eSLois Curfman McInnes       }
109204ca555eSLois Curfman McInnes       MPI_Allreduce((void*)&sum,(void*)norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);
109304ca555eSLois Curfman McInnes       *norm = sqrt(*norm);
109404ca555eSLois Curfman McInnes     }
109504ca555eSLois Curfman McInnes     else if (type == NORM_1) { /* max column norm */
109604ca555eSLois Curfman McInnes       double *tmp, *tmp2;
109704ca555eSLois Curfman McInnes       int    *jj, *garray = aij->garray;
109804ca555eSLois Curfman McInnes       tmp = (double *) PETSCMALLOC( aij->N*sizeof(double) ); CHKPTRQ(tmp);
109904ca555eSLois Curfman McInnes       tmp2 = (double *) PETSCMALLOC( aij->N*sizeof(double) ); CHKPTRQ(tmp2);
110004ca555eSLois Curfman McInnes       PETSCMEMSET(tmp,0,aij->N*sizeof(double));
110104ca555eSLois Curfman McInnes       *norm = 0.0;
110204ca555eSLois Curfman McInnes       v = amat->a; jj = amat->j;
110304ca555eSLois Curfman McInnes       for ( j=0; j<amat->nz; j++ ) {
110404ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
11056d531a12SLois Curfman McInnes         tmp[cstart + *jj++ - 1] += abs(*v++);
110604ca555eSLois Curfman McInnes #else
11076d531a12SLois Curfman McInnes         tmp[cstart + *jj++ - 1] += fabs(*v++);
110804ca555eSLois Curfman McInnes #endif
110904ca555eSLois Curfman McInnes       }
111004ca555eSLois Curfman McInnes       v = bmat->a; jj = bmat->j;
111104ca555eSLois Curfman McInnes       for ( j=0; j<bmat->nz; j++ ) {
111204ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
111304ca555eSLois Curfman McInnes         tmp[garray[*jj++ - 1]] += abs(*v++);
111404ca555eSLois Curfman McInnes #else
111504ca555eSLois Curfman McInnes         tmp[garray[*jj++ - 1]] += fabs(*v++);
111604ca555eSLois Curfman McInnes #endif
111704ca555eSLois Curfman McInnes       }
111804ca555eSLois Curfman McInnes       MPI_Allreduce((void*)tmp,(void*)tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);
111904ca555eSLois Curfman McInnes       for ( j=0; j<aij->N; j++ ) {
112004ca555eSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
112104ca555eSLois Curfman McInnes       }
112204ca555eSLois Curfman McInnes       PETSCFREE(tmp); PETSCFREE(tmp2);
112304ca555eSLois Curfman McInnes     }
112404ca555eSLois Curfman McInnes     else if (type == NORM_INFINITY) { /* max row norm */
112504ca555eSLois Curfman McInnes       double normtemp = 0.0;
112604ca555eSLois Curfman McInnes       for ( j=0; j<amat->m; j++ ) {
112704ca555eSLois Curfman McInnes         v = amat->a + amat->i[j] - 1;
112804ca555eSLois Curfman McInnes         sum = 0.0;
112904ca555eSLois Curfman McInnes         for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) {
113004ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
113104ca555eSLois Curfman McInnes           sum += abs(*v); v++;
113204ca555eSLois Curfman McInnes #else
113304ca555eSLois Curfman McInnes           sum += fabs(*v); v++;
113404ca555eSLois Curfman McInnes #endif
113504ca555eSLois Curfman McInnes         }
113604ca555eSLois Curfman McInnes         v = bmat->a + bmat->i[j] - 1;
113704ca555eSLois Curfman McInnes         for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) {
113804ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
113904ca555eSLois Curfman McInnes           sum += abs(*v); v++;
114004ca555eSLois Curfman McInnes #else
114104ca555eSLois Curfman McInnes           sum += fabs(*v); v++;
114204ca555eSLois Curfman McInnes #endif
114304ca555eSLois Curfman McInnes         }
114404ca555eSLois Curfman McInnes         if (sum > normtemp) normtemp = sum;
114504ca555eSLois Curfman McInnes         MPI_Allreduce((void*)&normtemp,(void*)norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);
114604ca555eSLois Curfman McInnes       }
114704ca555eSLois Curfman McInnes     }
114804ca555eSLois Curfman McInnes     else {
114904ca555eSLois Curfman McInnes       SETERRQ(1,"MatNorm_MPIRow:No support for the two norm");
115004ca555eSLois Curfman McInnes     }
115137fa93a5SLois Curfman McInnes   }
1152855ac2c5SLois Curfman McInnes   return 0;
1153855ac2c5SLois Curfman McInnes }
1154855ac2c5SLois Curfman McInnes 
11550de55854SLois Curfman McInnes static int MatTranspose_MPIAIJ(Mat A,Mat *matout)
1156b7c46309SBarry Smith {
1157b7c46309SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1158b7c46309SBarry Smith   int        ierr;
1159b7c46309SBarry Smith   Mat        B;
1160b7c46309SBarry Smith   Mat_AIJ    *Aloc;
1161b7c46309SBarry Smith   int        M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct;
1162b7c46309SBarry Smith   Scalar     *array;
1163b7c46309SBarry Smith 
11640de55854SLois Curfman McInnes   if (!matout && M != N) SETERRQ(1,
11650de55854SLois Curfman McInnes     "MatTranspose_MPIAIJ: Cannot transpose rectangular matrix in place");
1166b7c46309SBarry Smith   ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,0,0,0,&B);
1167b7c46309SBarry Smith   CHKERRQ(ierr);
1168b7c46309SBarry Smith 
1169b7c46309SBarry Smith   /* copy over the A part */
1170b7c46309SBarry Smith   Aloc = (Mat_AIJ*) a->A->data;
1171b7c46309SBarry Smith   m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1172b7c46309SBarry Smith   row = a->rstart;
1173b7c46309SBarry Smith   for ( i=0; i<ai[m]-1; i++ ) {aj[i] += a->cstart - 1;}
1174b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1175b7c46309SBarry Smith       ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERTVALUES);
1176b7c46309SBarry Smith       CHKERRQ(ierr);
1177b7c46309SBarry Smith       row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
1178b7c46309SBarry Smith   }
1179b7c46309SBarry Smith   aj = Aloc->j;
1180b7c46309SBarry Smith   for ( i=0; i<ai[m]-1; i++ ) {aj[i] -= a->cstart - 1;}
1181b7c46309SBarry Smith 
1182b7c46309SBarry Smith   /* copy over the B part */
1183b7c46309SBarry Smith   Aloc = (Mat_AIJ*) a->B->data;
1184b7c46309SBarry Smith   m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1185b7c46309SBarry Smith   row = a->rstart;
1186b7c46309SBarry Smith   ct = cols = (int *) PETSCMALLOC( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
1187b7c46309SBarry Smith   for ( i=0; i<ai[m]-1; i++ ) {cols[i] = a->garray[aj[i]-1];}
1188b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1189b7c46309SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERTVALUES);
1190b7c46309SBarry Smith     CHKERRQ(ierr);
1191b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
1192b7c46309SBarry Smith   }
1193b7c46309SBarry Smith   PETSCFREE(ct);
1194b7c46309SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1195b7c46309SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
11960de55854SLois Curfman McInnes   if (matout) {
11970de55854SLois Curfman McInnes     *matout = B;
11980de55854SLois Curfman McInnes   } else {
11990de55854SLois Curfman McInnes     /* This isn't really an in-place transpose .... but free data structures from a */
12000de55854SLois Curfman McInnes     PETSCFREE(a->rowners);
12010de55854SLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
12020de55854SLois Curfman McInnes     ierr = MatDestroy(a->B); CHKERRQ(ierr);
12030de55854SLois Curfman McInnes     if (a->colmap) PETSCFREE(a->colmap);
12040de55854SLois Curfman McInnes     if (a->garray) PETSCFREE(a->garray);
12050de55854SLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
12060de55854SLois Curfman McInnes     if (a->Mvctx) VecScatterCtxDestroy(a->Mvctx);
12070de55854SLois Curfman McInnes     PETSCFREE(a);
12080de55854SLois Curfman McInnes     PETSCMEMCPY(A,B,sizeof(struct _Mat));
12090de55854SLois Curfman McInnes     PETSCHEADERDESTROY(B);
12100de55854SLois Curfman McInnes   }
1211b7c46309SBarry Smith   return 0;
1212b7c46309SBarry Smith }
1213b7c46309SBarry Smith 
1214fc76ce87SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MatType,Mat *);
1215ff7509f2SBarry Smith static int MatCopyPrivate_MPIAIJ(Mat,Mat *);
1216d6dfbf8fSBarry Smith 
12178a729477SBarry Smith /* -------------------------------------------------------------------*/
12182ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
121939e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
122044a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
122144a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
1222299609e3SLois Curfman McInnes        MatSolve_MPIAIJ,MatSolveAdd_MPIAIJ,
1223299609e3SLois Curfman McInnes        MatSolveTrans_MPIAIJ,MatSolveTransAdd_MPIAIJ,
1224299609e3SLois Curfman McInnes        MatLUFactor_MPIAIJ,0,
122544a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
1226b7c46309SBarry Smith        MatTranspose_MPIAIJ,
1227a66be287SLois Curfman McInnes        MatGetInfo_MPIAIJ,0,
1228855ac2c5SLois Curfman McInnes        MatGetDiagonal_MPIAIJ,0,MatNorm_MPIAIJ,
1229ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
12301eb62cbbSBarry Smith        0,
1231299609e3SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,
1232299609e3SLois Curfman McInnes        MatGetReordering_MPIAIJ,
1233299609e3SLois Curfman McInnes        MatLUFactorSymbolic_MPIAIJ,MatLUFactorNumeric_MPIAIJ,0,0,
1234d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
1235299609e3SLois Curfman McInnes        MatILUFactorSymbolic_MPIAIJ,0,
1236ff7509f2SBarry Smith        0,0,MatConvert_MPIAIJ,0,0,MatCopyPrivate_MPIAIJ};
12378a729477SBarry Smith 
12388a729477SBarry Smith /*@
1239ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
1240ff756334SLois Curfman McInnes    (the default parallel PETSc format).
12418a729477SBarry Smith 
12428a729477SBarry Smith    Input Parameters:
12431eb62cbbSBarry Smith .  comm - MPI communicator
12447d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
12457d3e4905SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated
12467d3e4905SLois Curfman McInnes            if N is given)
12477d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
12487d3e4905SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated
12497d3e4905SLois Curfman McInnes            if n is given)
1250ab693e5aSLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of local submatrix
1251ff756334SLois Curfman McInnes            (same for all local rows)
1252ab693e5aSLois Curfman McInnes .  d_nzz - number of nonzeros per row in diagonal portion of local submatrix
1253ab693e5aSLois Curfman McInnes            or null (possibly different for each row).  You must leave room
1254ab693e5aSLois Curfman McInnes            for the diagonal entry even if it is zero.
1255ab693e5aSLois Curfman McInnes .  o_nz - number of nonzeros per row in off-diagonal portion of local
1256ab693e5aSLois Curfman McInnes            submatrix (same for all local rows).
1257ab693e5aSLois Curfman McInnes .  o_nzz - number of nonzeros per row in off-diagonal portion of local
1258ab693e5aSLois Curfman McInnes            submatrix or null (possibly different for each row).
12598a729477SBarry Smith 
1260ff756334SLois Curfman McInnes    Output Parameter:
12618a729477SBarry Smith .  newmat - the matrix
12628a729477SBarry Smith 
1263ff756334SLois Curfman McInnes    Notes:
1264ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1265ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
1266ff756334SLois Curfman McInnes    storage.  That is, the stored row and column indices begin at
1267ab693e5aSLois Curfman McInnes    one, not zero.  See the users manual for further details.
1268ff756334SLois Curfman McInnes 
1269ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1270ff756334SLois Curfman McInnes    (possibly both).
1271ff756334SLois Curfman McInnes 
1272e0245417SLois Curfman McInnes    Storage Information:
1273e0245417SLois Curfman McInnes    For a square global matrix we define each processor's diagonal portion
1274e0245417SLois Curfman McInnes    to be its local rows and the corresponding columns (a square submatrix);
1275e0245417SLois Curfman McInnes    each processor's off-diagonal portion encompasses the remainder of the
1276e0245417SLois Curfman McInnes    local matrix (a rectangular submatrix).
1277e0245417SLois Curfman McInnes 
1278e0245417SLois Curfman McInnes    The user can specify preallocated storage for the diagonal part of
1279e0245417SLois Curfman McInnes    the local submatrix with either d_nz or d_nnz (not both). Set both
1280e0245417SLois Curfman McInnes    d_nz and d_nnz to zero for PETSc to control dynamic memory allocation.
1281e0245417SLois Curfman McInnes    Likewise, specify preallocated storage for the off-diagonal part of
1282e0245417SLois Curfman McInnes    the local submatrix with o_nz or o_nnz (not both).
1283ff756334SLois Curfman McInnes 
1284dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel
1285ff756334SLois Curfman McInnes 
1286*fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
12878a729477SBarry Smith @*/
12881eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,
12891eb62cbbSBarry Smith                  int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat)
12908a729477SBarry Smith {
12918a729477SBarry Smith   Mat          mat;
129244a69424SLois Curfman McInnes   Mat_MPIAIJ   *aij;
12936abc6512SBarry Smith   int          ierr, i,sum[2],work[2];
12948a729477SBarry Smith   *newmat         = 0;
129544a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm);
1296a5a9c739SBarry Smith   PLogObjectCreate(mat);
129778b31e54SBarry Smith   mat->data       = (void *) (aij = PETSCNEW(Mat_MPIAIJ)); CHKPTRQ(aij);
12988a729477SBarry Smith   mat->ops        = &MatOps;
129944a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
130044a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
13018a729477SBarry Smith   mat->factor     = 0;
1302d6dfbf8fSBarry Smith 
130307a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
13041eb62cbbSBarry Smith   MPI_Comm_rank(comm,&aij->mytid);
13051eb62cbbSBarry Smith   MPI_Comm_size(comm,&aij->numtids);
13061eb62cbbSBarry Smith 
1307dbd7a890SLois Curfman McInnes   if (M == PETSC_DECIDE || N == PETSC_DECIDE) {
13081eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1309d6dfbf8fSBarry Smith     MPI_Allreduce((void *) work,(void *) sum,2,MPI_INT,MPI_SUM,comm );
1310dbd7a890SLois Curfman McInnes     if (M == PETSC_DECIDE) M = sum[0];
1311dbd7a890SLois Curfman McInnes     if (N == PETSC_DECIDE) N = sum[1];
13121eb62cbbSBarry Smith   }
1313dbd7a890SLois Curfman McInnes   if (m == PETSC_DECIDE)
1314dbd7a890SLois Curfman McInnes     {m = M/aij->numtids + ((M % aij->numtids) > aij->mytid);}
1315dbd7a890SLois Curfman McInnes   if (n == PETSC_DECIDE)
1316dbd7a890SLois Curfman McInnes     {n = N/aij->numtids + ((N % aij->numtids) > aij->mytid);}
13178a729477SBarry Smith   aij->m = m;
13188a729477SBarry Smith   aij->n = n;
13191eb62cbbSBarry Smith   aij->N = N;
13201eb62cbbSBarry Smith   aij->M = M;
13211eb62cbbSBarry Smith 
13221eb62cbbSBarry Smith   /* build local table of row and column ownerships */
132378b31e54SBarry Smith   aij->rowners = (int *) PETSCMALLOC(2*(aij->numtids+2)*sizeof(int));
132478b31e54SBarry Smith   CHKPTRQ(aij->rowners);
1325464493b3SBarry Smith   PLogObjectMemory(mat,2*(aij->numtids+2)*sizeof(int)+sizeof(struct _Mat)+
1326464493b3SBarry Smith                        sizeof(Mat_MPIAIJ));
13271eb62cbbSBarry Smith   aij->cowners = aij->rowners + aij->numtids + 1;
13281eb62cbbSBarry Smith   MPI_Allgather(&m,1,MPI_INT,aij->rowners+1,1,MPI_INT,comm);
13291eb62cbbSBarry Smith   aij->rowners[0] = 0;
13301eb62cbbSBarry Smith   for ( i=2; i<=aij->numtids; i++ ) {
13311eb62cbbSBarry Smith     aij->rowners[i] += aij->rowners[i-1];
13328a729477SBarry Smith   }
13331eb62cbbSBarry Smith   aij->rstart = aij->rowners[aij->mytid];
13341eb62cbbSBarry Smith   aij->rend   = aij->rowners[aij->mytid+1];
13351eb62cbbSBarry Smith   MPI_Allgather(&n,1,MPI_INT,aij->cowners+1,1,MPI_INT,comm);
13361eb62cbbSBarry Smith   aij->cowners[0] = 0;
13371eb62cbbSBarry Smith   for ( i=2; i<=aij->numtids; i++ ) {
13381eb62cbbSBarry Smith     aij->cowners[i] += aij->cowners[i-1];
13398a729477SBarry Smith   }
13401eb62cbbSBarry Smith   aij->cstart = aij->cowners[aij->mytid];
13411eb62cbbSBarry Smith   aij->cend   = aij->cowners[aij->mytid+1];
13428a729477SBarry Smith 
1343*fafbff53SBarry Smith   ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&aij->A);
134478b31e54SBarry Smith   CHKERRQ(ierr);
1345a5a9c739SBarry Smith   PLogObjectParent(mat,aij->A);
1346*fafbff53SBarry Smith   ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&aij->B);
134778b31e54SBarry Smith   CHKERRQ(ierr);
1348a5a9c739SBarry Smith   PLogObjectParent(mat,aij->B);
13498a729477SBarry Smith 
13501eb62cbbSBarry Smith   /* build cache for off array entries formed */
135178b31e54SBarry Smith   ierr = StashBuild_Private(&aij->stash); CHKERRQ(ierr);
13529e25ed09SBarry Smith   aij->colmap    = 0;
13539e25ed09SBarry Smith   aij->garray    = 0;
13548a729477SBarry Smith 
13551eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
13561eb62cbbSBarry Smith   aij->lvec      = 0;
13571eb62cbbSBarry Smith   aij->Mvctx     = 0;
1358d6dfbf8fSBarry Smith   aij->assembled = 0;
13598a729477SBarry Smith 
1360d6dfbf8fSBarry Smith   *newmat = mat;
1361d6dfbf8fSBarry Smith   return 0;
1362d6dfbf8fSBarry Smith }
1363c74985f6SBarry Smith 
1364ff7509f2SBarry Smith static int MatCopyPrivate_MPIAIJ(Mat matin,Mat *newmat)
1365d6dfbf8fSBarry Smith {
1366d6dfbf8fSBarry Smith   Mat        mat;
136744a69424SLois Curfman McInnes   Mat_MPIAIJ *aij,*oldmat = (Mat_MPIAIJ *) matin->data;
13682ee70a88SLois Curfman McInnes   int        ierr, len;
1369d6dfbf8fSBarry Smith   *newmat      = 0;
1370d6dfbf8fSBarry Smith 
1371bbb6d6a8SBarry Smith   if (!oldmat->assembled)
1372bbb6d6a8SBarry Smith     SETERRQ(1,"MatCopyPrivate_MPIAIJ:Cannot copy unassembled matrix");
137344a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm);
1374a5a9c739SBarry Smith   PLogObjectCreate(mat);
137578b31e54SBarry Smith   mat->data       = (void *) (aij = PETSCNEW(Mat_MPIAIJ)); CHKPTRQ(aij);
1376d6dfbf8fSBarry Smith   mat->ops        = &MatOps;
137744a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
137844a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
1379d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1380d6dfbf8fSBarry Smith 
1381d6dfbf8fSBarry Smith   aij->m          = oldmat->m;
1382d6dfbf8fSBarry Smith   aij->n          = oldmat->n;
1383d6dfbf8fSBarry Smith   aij->M          = oldmat->M;
1384d6dfbf8fSBarry Smith   aij->N          = oldmat->N;
1385d6dfbf8fSBarry Smith 
1386d6dfbf8fSBarry Smith   aij->assembled  = 1;
1387d6dfbf8fSBarry Smith   aij->rstart     = oldmat->rstart;
1388d6dfbf8fSBarry Smith   aij->rend       = oldmat->rend;
1389d6dfbf8fSBarry Smith   aij->cstart     = oldmat->cstart;
1390d6dfbf8fSBarry Smith   aij->cend       = oldmat->cend;
1391d6dfbf8fSBarry Smith   aij->numtids    = oldmat->numtids;
1392d6dfbf8fSBarry Smith   aij->mytid      = oldmat->mytid;
139307a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
1394d6dfbf8fSBarry Smith 
139578b31e54SBarry Smith   aij->rowners        = (int *) PETSCMALLOC( (aij->numtids+1)*sizeof(int) );
139678b31e54SBarry Smith   CHKPTRQ(aij->rowners);
1397464493b3SBarry Smith   PLogObjectMemory(mat,(aij->numtids+1)*sizeof(int)+sizeof(struct _Mat)+
1398464493b3SBarry Smith                        sizeof(Mat_MPIAIJ));
139978b31e54SBarry Smith   PETSCMEMCPY(aij->rowners,oldmat->rowners,(aij->numtids+1)*sizeof(int));
140078b31e54SBarry Smith   ierr = StashInitialize_Private(&aij->stash); CHKERRQ(ierr);
14012ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
140278b31e54SBarry Smith     aij->colmap      = (int *) PETSCMALLOC( (aij->N)*sizeof(int) );
140378b31e54SBarry Smith     CHKPTRQ(aij->colmap);
1404464493b3SBarry Smith     PLogObjectMemory(mat,(aij->N)*sizeof(int));
140578b31e54SBarry Smith     PETSCMEMCPY(aij->colmap,oldmat->colmap,(aij->N)*sizeof(int));
14062ee70a88SLois Curfman McInnes   } else aij->colmap = 0;
14072ee70a88SLois Curfman McInnes   if (oldmat->garray && (len = ((Mat_AIJ *) (oldmat->B->data))->n)) {
140878b31e54SBarry Smith     aij->garray      = (int *) PETSCMALLOC(len*sizeof(int) ); CHKPTRQ(aij->garray);
1409464493b3SBarry Smith     PLogObjectMemory(mat,len*sizeof(int));
141078b31e54SBarry Smith     PETSCMEMCPY(aij->garray,oldmat->garray,len*sizeof(int));
14112ee70a88SLois Curfman McInnes   } else aij->garray = 0;
1412d6dfbf8fSBarry Smith 
141378b31e54SBarry Smith   ierr =  VecDuplicate(oldmat->lvec,&aij->lvec); CHKERRQ(ierr);
1414a5a9c739SBarry Smith   PLogObjectParent(mat,aij->lvec);
141578b31e54SBarry Smith   ierr =  VecScatterCtxCopy(oldmat->Mvctx,&aij->Mvctx); CHKERRQ(ierr);
1416a5a9c739SBarry Smith   PLogObjectParent(mat,aij->Mvctx);
141778b31e54SBarry Smith   ierr =  MatConvert(oldmat->A,MATSAME,&aij->A); CHKERRQ(ierr);
1418a5a9c739SBarry Smith   PLogObjectParent(mat,aij->A);
141978b31e54SBarry Smith   ierr =  MatConvert(oldmat->B,MATSAME,&aij->B); CHKERRQ(ierr);
1420a5a9c739SBarry Smith   PLogObjectParent(mat,aij->B);
14218a729477SBarry Smith   *newmat = mat;
14228a729477SBarry Smith   return 0;
14238a729477SBarry Smith }
1424