xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision 80b4ade8ca6f14c7a8695353deb39be397bcc1ce)
1cb512458SBarry Smith #ifndef lint
2*80b4ade8SLois Curfman McInnes static char vcid[] = "$Id: mpiaij.c,v 1.63 1995/08/15 20:28:20 bsmith Exp curfman $";
3cb512458SBarry Smith #endif
48a729477SBarry Smith 
51eb62cbbSBarry Smith #include "mpiaij.h"
68a729477SBarry Smith #include "vec/vecimpl.h"
7d6dfbf8fSBarry Smith #include "inline/spops.h"
88a729477SBarry Smith 
99e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column
109e25ed09SBarry Smith number to the local number in the off-diagonal part of the local
119e25ed09SBarry Smith storage of the matrix.  This is done in a non scable way since the
129e25ed09SBarry Smith length of colmap equals the global matrix length.
139e25ed09SBarry Smith */
14b7c46309SBarry Smith static int CreateColmap_Private(Mat mat)
159e25ed09SBarry Smith {
1644a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
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 */
3916b5873e3SBarry Smith   ierr = ISCreateSequential(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)
500a5a9c739SBarry Smith   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"
515ee50ffe9SBarry Smith #include "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;
5211eb62cbbSBarry Smith   int        ierr;
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 */
526154123eaSLois Curfman McInnes     viewer = STDOUT_VIEWER; vobj = (PetscObject) viewer;
527154123eaSLois Curfman McInnes   }
528ab254492SBarry Smith   if (vobj->cookie == DRAW_COOKIE && vobj->type == NULLWINDOW) return 0;
5297857610eSBarry Smith   if (vobj->cookie == VIEWER_COOKIE && vobj->type == FILE_VIEWER) {
5304a0ce102SLois Curfman McInnes     FILE *fd = ViewerFileGetPointer_Private(viewer);
5317f813858SBarry Smith     MPIU_Seq_begin(mat->comm,1);
532d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
5331eb62cbbSBarry Smith              aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
5341eb62cbbSBarry Smith              aij->cend);
53578b31e54SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
53678b31e54SBarry Smith     ierr = MatView(aij->B,viewer); CHKERRQ(ierr);
537d13ab20cSBarry Smith     fflush(fd);
5387f813858SBarry Smith     MPIU_Seq_end(mat->comm,1);
539d13ab20cSBarry Smith   }
5407857610eSBarry Smith   else if ((vobj->cookie == VIEWER_COOKIE && vobj->type == FILES_VIEWER) ||
5417857610eSBarry Smith             vobj->cookie == DRAW_COOKIE) {
54295373324SBarry Smith     int numtids = aij->numtids, mytid = aij->mytid;
54395373324SBarry Smith     if (numtids == 1) {
54478b31e54SBarry Smith       ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
54595373324SBarry Smith     }
54695373324SBarry Smith     else {
54795373324SBarry Smith       /* assemble the entire matrix onto first processor. */
54895373324SBarry Smith       Mat     A;
5492ee70a88SLois Curfman McInnes       Mat_AIJ *Aloc;
5502eb8c8abSBarry Smith       int     M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
55195373324SBarry Smith       Scalar  *a;
5522ee70a88SLois Curfman McInnes 
55395373324SBarry Smith       if (!mytid) {
55495373324SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,0,0,0,&A);
55595373324SBarry Smith       }
55695373324SBarry Smith       else {
55795373324SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,0,0,0,&A);
55895373324SBarry Smith       }
559464493b3SBarry Smith       PLogObjectParent(mat,A);
56078b31e54SBarry Smith       CHKERRQ(ierr);
56195373324SBarry Smith 
56295373324SBarry Smith       /* copy over the A part */
5632ee70a88SLois Curfman McInnes       Aloc = (Mat_AIJ*) aij->A->data;
5642ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
56595373324SBarry Smith       row = aij->rstart;
566b7c46309SBarry Smith       for ( i=0; i<ai[m]-1; i++ ) {aj[i] += aij->cstart - 1;}
56795373324SBarry Smith       for ( i=0; i<m; i++ ) {
56807a0e7adSLois Curfman McInnes         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERTVALUES);
56978b31e54SBarry Smith         CHKERRQ(ierr);
57095373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
57195373324SBarry Smith       }
5722ee70a88SLois Curfman McInnes       aj = Aloc->j;
573b7c46309SBarry Smith       for ( i=0; i<ai[m]-1; i++ ) {aj[i] -= aij->cstart - 1;}
57495373324SBarry Smith 
57595373324SBarry Smith       /* copy over the B part */
5762ee70a88SLois Curfman McInnes       Aloc = (Mat_AIJ*) aij->B->data;
5772ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
57895373324SBarry Smith       row = aij->rstart;
57978b31e54SBarry Smith       ct = cols = (int *) PETSCMALLOC( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
580b7c46309SBarry Smith       for ( i=0; i<ai[m]-1; i++ ) {cols[i] = aij->garray[aj[i]-1];}
58195373324SBarry Smith       for ( i=0; i<m; i++ ) {
58207a0e7adSLois Curfman McInnes         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERTVALUES);
58378b31e54SBarry Smith         CHKERRQ(ierr);
58495373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
58595373324SBarry Smith       }
58678b31e54SBarry Smith       PETSCFREE(ct);
58778b31e54SBarry Smith       ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
58878b31e54SBarry Smith       ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
58995373324SBarry Smith       if (!mytid) {
59078b31e54SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr);
59195373324SBarry Smith       }
59278b31e54SBarry Smith       ierr = MatDestroy(A); CHKERRQ(ierr);
59395373324SBarry Smith     }
59495373324SBarry Smith   }
5951eb62cbbSBarry Smith   return 0;
5961eb62cbbSBarry Smith }
5971eb62cbbSBarry Smith 
598*80b4ade8SLois Curfman McInnes extern int MatMarkDiag_AIJ(Mat);
5991eb62cbbSBarry Smith /*
6001eb62cbbSBarry Smith     This has to provide several versions.
6011eb62cbbSBarry Smith 
6021eb62cbbSBarry Smith      1) per sequential
6031eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
6041eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
605d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
6061eb62cbbSBarry Smith */
607fc76ce87SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
608fc76ce87SLois Curfman McInnes                            double shift,int its,Vec xx)
6098a729477SBarry Smith {
61044a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
611d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
61244a69424SLois Curfman McInnes   Mat_AIJ    *A = (Mat_AIJ *) AA->data, *B = (Mat_AIJ *)BB->data;
6136abc6512SBarry Smith   Scalar     zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts;
6146abc6512SBarry Smith   int        ierr,*idx, *diag;
6156abc6512SBarry Smith   int        n = mat->n, m = mat->m, i;
616da3a660dSBarry Smith   Vec        tt;
6178a729477SBarry Smith 
61878b31e54SBarry Smith   if (!mat->assembled) SETERRQ(1,"MatRelax_MPIAIJ:must assemble matrix first");
6191eb62cbbSBarry Smith 
620d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
621d6dfbf8fSBarry Smith   xs = x -1; /* shift by one for index start of 1 */
622da3a660dSBarry Smith   ls--;
623*80b4ade8SLois Curfman McInnes   if (!A->diag) {if ((ierr = MatMarkDiag_AIJ(AA))) return ierr;}
624d6dfbf8fSBarry Smith   diag = A->diag;
625acb40c82SBarry Smith   if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) {
62642f8db3fSLois Curfman McInnes     SETERRQ(1,"MatRelax_MPIAIJ:Option not yet supported");
627acb40c82SBarry Smith   }
628da3a660dSBarry Smith   if (flag & SOR_EISENSTAT) {
629da3a660dSBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
630da3a660dSBarry Smith     U is upper triangular, E is diagonal; This routine applies
631da3a660dSBarry Smith 
632da3a660dSBarry Smith             (L + E)^{-1} A (U + E)^{-1}
633da3a660dSBarry Smith 
634da3a660dSBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
635da3a660dSBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
636da3a660dSBarry Smith     is the relaxation factor.
637da3a660dSBarry Smith     */
63878b31e54SBarry Smith     ierr = VecDuplicate(xx,&tt); CHKERRQ(ierr);
639464493b3SBarry Smith     PLogObjectParent(matin,tt);
640da3a660dSBarry Smith     VecGetArray(tt,&t);
641da3a660dSBarry Smith     scale = (2.0/omega) - 1.0;
642da3a660dSBarry Smith     /*  x = (E + U)^{-1} b */
643da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
64406be10caSBarry Smith     ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
64578b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
646da3a660dSBarry Smith     for ( i=m-1; i>-1; i-- ) {
647da3a660dSBarry Smith       n    = A->i[i+1] - diag[i] - 1;
648da3a660dSBarry Smith       idx  = A->j + diag[i];
649da3a660dSBarry Smith       v    = A->a + diag[i];
650da3a660dSBarry Smith       sum  = b[i];
651da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
652da3a660dSBarry Smith       d    = shift + A->a[diag[i]-1];
653da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
654da3a660dSBarry Smith       idx  = B->j + B->i[i] - 1;
655da3a660dSBarry Smith       v    = B->a + B->i[i] - 1;
656da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
657da3a660dSBarry Smith       x[i] = omega*(sum/d);
658da3a660dSBarry Smith     }
659edd2f0e1SBarry Smith     ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
66078b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
661da3a660dSBarry Smith 
662da3a660dSBarry Smith     /*  t = b - (2*E - D)x */
663da3a660dSBarry Smith     v = A->a;
664da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ - 1])*x[i]; }
665da3a660dSBarry Smith 
666da3a660dSBarry Smith     /*  t = (E + L)^{-1}t */
667da3a660dSBarry Smith     ts = t - 1; /* shifted by one for index start of a or mat->j*/
668da3a660dSBarry Smith     diag = A->diag;
669da3a660dSBarry Smith     VecSet(&zero,mat->lvec);
67006be10caSBarry Smith     ierr = VecPipelineBegin(tt,mat->lvec,INSERTVALUES,PIPELINEDOWN,
67178b31e54SBarry Smith                                                  mat->Mvctx); CHKERRQ(ierr);
672da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
673da3a660dSBarry Smith       n    = diag[i] - A->i[i];
674da3a660dSBarry Smith       idx  = A->j + A->i[i] - 1;
675da3a660dSBarry Smith       v    = A->a + A->i[i] - 1;
676da3a660dSBarry Smith       sum  = t[i];
677da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
678da3a660dSBarry Smith       d    = shift + A->a[diag[i]-1];
679da3a660dSBarry Smith       n    = B->i[i+1] - B->i[i];
680da3a660dSBarry Smith       idx  = B->j + B->i[i] - 1;
681da3a660dSBarry Smith       v    = B->a + B->i[i] - 1;
682da3a660dSBarry Smith       SPARSEDENSEMDOT(sum,ls,v,idx,n);
683da3a660dSBarry Smith       t[i] = omega*(sum/d);
684da3a660dSBarry Smith     }
685edd2f0e1SBarry Smith     ierr = VecPipelineEnd(tt,mat->lvec,INSERTVALUES,PIPELINEDOWN,
68678b31e54SBarry Smith                                                     mat->Mvctx); CHKERRQ(ierr);
687da3a660dSBarry Smith     /*  x = x + t */
688da3a660dSBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
689da3a660dSBarry Smith     VecDestroy(tt);
690da3a660dSBarry Smith     return 0;
691da3a660dSBarry Smith   }
692da3a660dSBarry Smith 
6931eb62cbbSBarry Smith 
694d6dfbf8fSBarry Smith   if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){
695da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
696da3a660dSBarry Smith       VecSet(&zero,mat->lvec); VecSet(&zero,xx);
697da3a660dSBarry Smith     }
698da3a660dSBarry Smith     else {
69906be10caSBarry Smith       ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
70078b31e54SBarry Smith       CHKERRQ(ierr);
70106be10caSBarry Smith       ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
70278b31e54SBarry Smith       CHKERRQ(ierr);
703da3a660dSBarry Smith     }
704d6dfbf8fSBarry Smith     while (its--) {
705d6dfbf8fSBarry Smith       /* go down through the rows */
70606be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
70778b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
708d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
709d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
710d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
711d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
712d6dfbf8fSBarry Smith         sum  = b[i];
713d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
714d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
715d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
716d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
717d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
718d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
719d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
720d6dfbf8fSBarry Smith       }
721edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
72278b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
723d6dfbf8fSBarry Smith       /* come up through the rows */
72406be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
72578b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
726d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
727d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
728d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
729d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
730d6dfbf8fSBarry Smith         sum  = b[i];
731d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
732d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
733d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
734d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
735d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
736d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
737d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
738d6dfbf8fSBarry Smith       }
739edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
74078b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
741d6dfbf8fSBarry Smith     }
742d6dfbf8fSBarry Smith   }
743d6dfbf8fSBarry Smith   else if (flag & SOR_FORWARD_SWEEP){
744da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
745da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
74606be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
74778b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
748da3a660dSBarry Smith       for ( i=0; i<m; i++ ) {
749da3a660dSBarry Smith         n    = diag[i] - A->i[i];
750da3a660dSBarry Smith         idx  = A->j + A->i[i] - 1;
751da3a660dSBarry Smith         v    = A->a + A->i[i] - 1;
752da3a660dSBarry Smith         sum  = b[i];
753da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
754da3a660dSBarry Smith         d    = shift + A->a[diag[i]-1];
755da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
756da3a660dSBarry Smith         idx  = B->j + B->i[i] - 1;
757da3a660dSBarry Smith         v    = B->a + B->i[i] - 1;
758da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
759da3a660dSBarry Smith         x[i] = omega*(sum/d);
760da3a660dSBarry Smith       }
761edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
76278b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
763da3a660dSBarry Smith       its--;
764da3a660dSBarry Smith     }
765d6dfbf8fSBarry Smith     while (its--) {
76606be10caSBarry Smith       ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
76778b31e54SBarry Smith       CHKERRQ(ierr);
76806be10caSBarry Smith       ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERUP,mat->Mvctx);
76978b31e54SBarry Smith       CHKERRQ(ierr);
77006be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
77178b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
772d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
773d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
774d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
775d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
776d6dfbf8fSBarry Smith         sum  = b[i];
777d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
778d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
779d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
780d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
781d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
782d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
783d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
784d6dfbf8fSBarry Smith       }
785edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEDOWN,
78678b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
787d6dfbf8fSBarry Smith     }
788d6dfbf8fSBarry Smith   }
789d6dfbf8fSBarry Smith   else if (flag & SOR_BACKWARD_SWEEP){
790da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
791da3a660dSBarry Smith       VecSet(&zero,mat->lvec);
79206be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
79378b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
794da3a660dSBarry Smith       for ( i=m-1; i>-1; i-- ) {
795da3a660dSBarry Smith         n    = A->i[i+1] - diag[i] - 1;
796da3a660dSBarry Smith         idx  = A->j + diag[i];
797da3a660dSBarry Smith         v    = A->a + diag[i];
798da3a660dSBarry Smith         sum  = b[i];
799da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
800da3a660dSBarry Smith         d    = shift + A->a[diag[i]-1];
801da3a660dSBarry Smith         n    = B->i[i+1] - B->i[i];
802da3a660dSBarry Smith         idx  = B->j + B->i[i] - 1;
803da3a660dSBarry Smith         v    = B->a + B->i[i] - 1;
804da3a660dSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
805da3a660dSBarry Smith         x[i] = omega*(sum/d);
806da3a660dSBarry Smith       }
807edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
80878b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
809da3a660dSBarry Smith       its--;
810da3a660dSBarry Smith     }
811d6dfbf8fSBarry Smith     while (its--) {
81206be10caSBarry Smith       ierr = VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERDOWN,
81378b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
81406be10caSBarry Smith       ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERDOWN,
81578b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
81606be10caSBarry Smith       ierr = VecPipelineBegin(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
81778b31e54SBarry Smith                               mat->Mvctx); CHKERRQ(ierr);
818d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
819d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
820d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
821d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
822d6dfbf8fSBarry Smith         sum  = b[i];
823d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
824d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
825d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
826d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
827d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
828d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
829d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
830d6dfbf8fSBarry Smith       }
831edd2f0e1SBarry Smith       ierr = VecPipelineEnd(xx,mat->lvec,INSERTVALUES,PIPELINEUP,
83278b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
833d6dfbf8fSBarry Smith     }
834d6dfbf8fSBarry Smith   }
835d6dfbf8fSBarry Smith   else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
836da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
837da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
838da3a660dSBarry Smith     }
83906be10caSBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
84078b31e54SBarry Smith     CHKERRQ(ierr);
84106be10caSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
84278b31e54SBarry Smith     CHKERRQ(ierr);
843d6dfbf8fSBarry Smith     while (its--) {
844d6dfbf8fSBarry Smith       /* go down through the rows */
845d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
846d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
847d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
848d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
849d6dfbf8fSBarry Smith         sum  = b[i];
850d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
851d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
852d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
853d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
854d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
855d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
856d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
857d6dfbf8fSBarry Smith       }
858d6dfbf8fSBarry Smith       /* come up through the rows */
859d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
860d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
861d6dfbf8fSBarry Smith         idx  = A->j + A->i[i] - 1;
862d6dfbf8fSBarry Smith         v    = A->a + A->i[i] - 1;
863d6dfbf8fSBarry Smith         sum  = b[i];
864d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
865d6dfbf8fSBarry Smith         d    = shift + A->a[diag[i]-1];
866d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
867d6dfbf8fSBarry Smith         idx  = B->j + B->i[i] - 1;
868d6dfbf8fSBarry Smith         v    = B->a + B->i[i] - 1;
869d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
870d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
871d6dfbf8fSBarry Smith       }
872d6dfbf8fSBarry Smith     }
873d6dfbf8fSBarry Smith   }
874d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_FORWARD_SWEEP){
875da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
876da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
877da3a660dSBarry Smith     }
87806be10caSBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
87978b31e54SBarry Smith     CHKERRQ(ierr);
88006be10caSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERALL,mat->Mvctx);
88178b31e54SBarry Smith     CHKERRQ(ierr);
882d6dfbf8fSBarry Smith     while (its--) {
883d6dfbf8fSBarry Smith       for ( i=0; i<m; 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_BACKWARD_SWEEP){
899da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
900da3a660dSBarry Smith       return MatRelax(mat->A,bb,omega,flag,shift,its,xx);
901da3a660dSBarry Smith     }
90206be10caSBarry Smith     ierr = VecScatterBegin(xx,mat->lvec,INSERTVALUES,SCATTERALL,
90378b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
90406be10caSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERTVALUES,SCATTERALL,
90578b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
906d6dfbf8fSBarry Smith     while (its--) {
907d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; 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   }
9228a729477SBarry Smith   return 0;
9238a729477SBarry Smith }
924a66be287SLois Curfman McInnes 
925fc76ce87SLois Curfman McInnes static int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,int *nz,
926fc76ce87SLois Curfman McInnes                              int *nzalloc,int *mem)
927a66be287SLois Curfman McInnes {
928a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
929a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
930a66be287SLois Curfman McInnes   int        ierr, isend[3], irecv[3], nzA, nzallocA, memA;
931a66be287SLois Curfman McInnes 
93278b31e54SBarry Smith   ierr = MatGetInfo(A,MAT_LOCAL,&nzA,&nzallocA,&memA); CHKERRQ(ierr);
93378b31e54SBarry Smith   ierr = MatGetInfo(B,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr);
934a66be287SLois Curfman McInnes   isend[0] += nzA; isend[1] += nzallocA; isend[2] += memA;
935a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
936a66be287SLois Curfman McInnes     *nz = isend[0]; *nzalloc = isend[1]; *mem = isend[2];
937a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
938a66be287SLois Curfman McInnes     MPI_Allreduce((void *) isend,(void *) irecv,3,MPI_INT,MPI_MAX,matin->comm);
939a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
940a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
941a66be287SLois Curfman McInnes     MPI_Allreduce((void *) isend,(void *) irecv,3,MPI_INT,MPI_SUM,matin->comm);
942a66be287SLois Curfman McInnes     *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2];
943a66be287SLois Curfman McInnes   }
944a66be287SLois Curfman McInnes   return 0;
945a66be287SLois Curfman McInnes }
946a66be287SLois Curfman McInnes 
947299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*);
948299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*);
949299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double);
950299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *);
951299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec);
952299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec);
953299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec);
954299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec);
955299609e3SLois Curfman McInnes 
956fc76ce87SLois Curfman McInnes static int MatSetOption_MPIAIJ(Mat aijin,MatOption op)
957c74985f6SBarry Smith {
95844a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data;
959c74985f6SBarry Smith 
960c74985f6SBarry Smith   if      (op == NO_NEW_NONZERO_LOCATIONS)  {
961c74985f6SBarry Smith     MatSetOption(aij->A,op);
962c74985f6SBarry Smith     MatSetOption(aij->B,op);
963c74985f6SBarry Smith   }
964c74985f6SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) {
965c74985f6SBarry Smith     MatSetOption(aij->A,op);
966c74985f6SBarry Smith     MatSetOption(aij->B,op);
967c74985f6SBarry Smith   }
968bbb6d6a8SBarry Smith   else if (op == COLUMN_ORIENTED)
969bbb6d6a8SBarry Smith     SETERRQ(1,"MatSetOption_MPIAIJ:Column oriented not supported");
970c74985f6SBarry Smith   return 0;
971c74985f6SBarry Smith }
972c74985f6SBarry Smith 
973d1710a03SLois Curfman McInnes static int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
974c74985f6SBarry Smith {
97544a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
976c74985f6SBarry Smith   *m = mat->M; *n = mat->N;
977c74985f6SBarry Smith   return 0;
978c74985f6SBarry Smith }
979c74985f6SBarry Smith 
980d1710a03SLois Curfman McInnes static int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
981c74985f6SBarry Smith {
98244a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
983b7c46309SBarry Smith   *m = mat->m; *n = mat->N;
984c74985f6SBarry Smith   return 0;
985c74985f6SBarry Smith }
986c74985f6SBarry Smith 
987d1710a03SLois Curfman McInnes static int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
988c74985f6SBarry Smith {
98944a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
990c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
991c74985f6SBarry Smith   return 0;
992c74985f6SBarry Smith }
993c74985f6SBarry Smith 
99439e00950SLois Curfman McInnes static int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
99539e00950SLois Curfman McInnes {
996154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
997154123eaSLois Curfman McInnes   Scalar     *vworkA, *vworkB, **pvA, **pvB;
998154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
999154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
100039e00950SLois Curfman McInnes 
100139e00950SLois Curfman McInnes   if (row < rstart || row >= rend)
1002bbb6d6a8SBarry Smith     SETERRQ(1,"MatGetRow_MPIAIJ:Currently you can get only local rows")
1003abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
100439e00950SLois Curfman McInnes 
1005154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1006154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1007154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
100878b31e54SBarry Smith   ierr = MatGetRow(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
100978b31e54SBarry Smith   ierr = MatGetRow(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
1010154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1011154123eaSLois Curfman McInnes 
1012154123eaSLois Curfman McInnes   if (v  || idx) {
1013154123eaSLois Curfman McInnes     if (nztot) {
1014154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
1015299609e3SLois Curfman McInnes       int imark;
101648b35521SBarry Smith       if (mat->assembled) {
1017154123eaSLois Curfman McInnes         for (i=0; i<nzB; i++) cworkB[i] = mat->garray[cworkB[i]];
101848b35521SBarry Smith       }
1019154123eaSLois Curfman McInnes       if (v) {
102078b31e54SBarry Smith         *v = (Scalar *) PETSCMALLOC( (nztot)*sizeof(Scalar) ); CHKPTRQ(*v);
102139e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1022154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*v)[i] = vworkB[i];
1023154123eaSLois Curfman McInnes           else break;
1024154123eaSLois Curfman McInnes         }
1025154123eaSLois Curfman McInnes         imark = i;
1026154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*v)[imark+i] = vworkA[i];
1027299609e3SLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*v)[nzA+i] = vworkB[i];
1028154123eaSLois Curfman McInnes       }
1029154123eaSLois Curfman McInnes       if (idx) {
103078b31e54SBarry Smith         *idx = (int *) PETSCMALLOC( (nztot)*sizeof(int) ); CHKPTRQ(*idx);
1031154123eaSLois Curfman McInnes         for (i=0; i<nzA; i++) cworkA[i] += cstart;
1032154123eaSLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
1033154123eaSLois Curfman McInnes           if (cworkB[i] < cstart)   (*idx)[i] = cworkB[i];
1034154123eaSLois Curfman McInnes           else break;
1035154123eaSLois Curfman McInnes         }
1036154123eaSLois Curfman McInnes         imark = i;
1037154123eaSLois Curfman McInnes         for ( i=0; i<nzA; i++ )     (*idx)[imark+i] = cworkA[i];
1038299609e3SLois Curfman McInnes         for ( i=imark; i<nzB; i++ ) (*idx)[nzA+i] = cworkB[i];
103939e00950SLois Curfman McInnes       }
104039e00950SLois Curfman McInnes     }
104139e00950SLois Curfman McInnes     else {*idx = 0; *v=0;}
1042154123eaSLois Curfman McInnes   }
104339e00950SLois Curfman McInnes   *nz = nztot;
104478b31e54SBarry Smith   ierr = MatRestoreRow(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
104578b31e54SBarry Smith   ierr = MatRestoreRow(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
104639e00950SLois Curfman McInnes   return 0;
104739e00950SLois Curfman McInnes }
104839e00950SLois Curfman McInnes 
104939e00950SLois Curfman McInnes static int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
105039e00950SLois Curfman McInnes {
105178b31e54SBarry Smith   if (idx) PETSCFREE(*idx);
105278b31e54SBarry Smith   if (v) PETSCFREE(*v);
105339e00950SLois Curfman McInnes   return 0;
105439e00950SLois Curfman McInnes }
105539e00950SLois Curfman McInnes 
1056855ac2c5SLois Curfman McInnes static int MatNorm_MPIAIJ(Mat mat,MatNormType type,double *norm)
1057855ac2c5SLois Curfman McInnes {
1058855ac2c5SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1059855ac2c5SLois Curfman McInnes   int ierr;
1060855ac2c5SLois Curfman McInnes   if (aij->numtids == 1) {
106114183eadSLois Curfman McInnes     ierr =  MatNorm(aij->A,type,norm); CHKERRQ(ierr);
1062855ac2c5SLois Curfman McInnes   } else
1063bbb6d6a8SBarry Smith     SETERRQ(1,"MatNorm_MPIAIJ:not yet supported in parallel");
1064855ac2c5SLois Curfman McInnes   return 0;
1065855ac2c5SLois Curfman McInnes }
1066855ac2c5SLois Curfman McInnes 
1067b7c46309SBarry Smith static int MatTranspose_MPIAIJ(Mat A,Mat *Bin)
1068b7c46309SBarry Smith {
1069b7c46309SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1070b7c46309SBarry Smith   int        ierr;
1071b7c46309SBarry Smith   Mat        B;
1072b7c46309SBarry Smith   Mat_AIJ    *Aloc;
1073b7c46309SBarry Smith   int        M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct;
1074b7c46309SBarry Smith   Scalar     *array;
1075b7c46309SBarry Smith 
1076b7c46309SBarry Smith   ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,0,0,0,&B);
1077b7c46309SBarry Smith   CHKERRQ(ierr);
1078b7c46309SBarry Smith 
1079b7c46309SBarry Smith   /* copy over the A part */
1080b7c46309SBarry Smith   Aloc = (Mat_AIJ*) a->A->data;
1081b7c46309SBarry Smith   m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1082b7c46309SBarry Smith   row = a->rstart;
1083b7c46309SBarry Smith   for ( i=0; i<ai[m]-1; i++ ) {aj[i] += a->cstart - 1;}
1084b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1085b7c46309SBarry Smith       ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERTVALUES);
1086b7c46309SBarry Smith       CHKERRQ(ierr);
1087b7c46309SBarry Smith       row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
1088b7c46309SBarry Smith   }
1089b7c46309SBarry Smith   aj = Aloc->j;
1090b7c46309SBarry Smith   for ( i=0; i<ai[m]-1; i++ ) {aj[i] -= a->cstart - 1;}
1091b7c46309SBarry Smith 
1092b7c46309SBarry Smith   /* copy over the B part */
1093b7c46309SBarry Smith   Aloc = (Mat_AIJ*) a->B->data;
1094b7c46309SBarry Smith   m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1095b7c46309SBarry Smith   row = a->rstart;
1096b7c46309SBarry Smith   ct = cols = (int *) PETSCMALLOC( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
1097b7c46309SBarry Smith   for ( i=0; i<ai[m]-1; i++ ) {cols[i] = a->garray[aj[i]-1];}
1098b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1099b7c46309SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERTVALUES);
1100b7c46309SBarry Smith     CHKERRQ(ierr);
1101b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
1102b7c46309SBarry Smith   }
1103b7c46309SBarry Smith   PETSCFREE(ct);
1104b7c46309SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1105b7c46309SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1106b7c46309SBarry Smith   *Bin = B;
1107b7c46309SBarry Smith   return 0;
1108b7c46309SBarry Smith }
1109b7c46309SBarry Smith 
1110fc76ce87SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MatType,Mat *);
1111ff7509f2SBarry Smith static int MatCopyPrivate_MPIAIJ(Mat,Mat *);
1112d6dfbf8fSBarry Smith 
11138a729477SBarry Smith /* -------------------------------------------------------------------*/
11142ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
111539e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
111644a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
111744a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
1118299609e3SLois Curfman McInnes        MatSolve_MPIAIJ,MatSolveAdd_MPIAIJ,
1119299609e3SLois Curfman McInnes        MatSolveTrans_MPIAIJ,MatSolveTransAdd_MPIAIJ,
1120299609e3SLois Curfman McInnes        MatLUFactor_MPIAIJ,0,
112144a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
1122b7c46309SBarry Smith        MatTranspose_MPIAIJ,
1123a66be287SLois Curfman McInnes        MatGetInfo_MPIAIJ,0,
1124855ac2c5SLois Curfman McInnes        MatGetDiagonal_MPIAIJ,0,MatNorm_MPIAIJ,
1125ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
11261eb62cbbSBarry Smith        0,
1127299609e3SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,
1128299609e3SLois Curfman McInnes        MatGetReordering_MPIAIJ,
1129299609e3SLois Curfman McInnes        MatLUFactorSymbolic_MPIAIJ,MatLUFactorNumeric_MPIAIJ,0,0,
1130d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
1131299609e3SLois Curfman McInnes        MatILUFactorSymbolic_MPIAIJ,0,
1132ff7509f2SBarry Smith        0,0,MatConvert_MPIAIJ,0,0,MatCopyPrivate_MPIAIJ};
11338a729477SBarry Smith 
11348a729477SBarry Smith /*@
1135ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
1136ff756334SLois Curfman McInnes    (the default parallel PETSc format).
11378a729477SBarry Smith 
11388a729477SBarry Smith    Input Parameters:
11391eb62cbbSBarry Smith .  comm - MPI communicator
11407d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
11417d3e4905SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated
11427d3e4905SLois Curfman McInnes            if N is given)
11437d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
11447d3e4905SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated
11457d3e4905SLois Curfman McInnes            if n is given)
1146ab693e5aSLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of local submatrix
1147ff756334SLois Curfman McInnes            (same for all local rows)
1148ab693e5aSLois Curfman McInnes .  d_nzz - number of nonzeros per row in diagonal portion of local submatrix
1149ab693e5aSLois Curfman McInnes            or null (possibly different for each row).  You must leave room
1150ab693e5aSLois Curfman McInnes            for the diagonal entry even if it is zero.
1151ab693e5aSLois Curfman McInnes .  o_nz - number of nonzeros per row in off-diagonal portion of local
1152ab693e5aSLois Curfman McInnes            submatrix (same for all local rows).
1153ab693e5aSLois Curfman McInnes .  o_nzz - number of nonzeros per row in off-diagonal portion of local
1154ab693e5aSLois Curfman McInnes            submatrix or null (possibly different for each row).
11558a729477SBarry Smith 
1156ff756334SLois Curfman McInnes    Output Parameter:
11578a729477SBarry Smith .  newmat - the matrix
11588a729477SBarry Smith 
1159ff756334SLois Curfman McInnes    Notes:
1160ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1161ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
1162ff756334SLois Curfman McInnes    storage.  That is, the stored row and column indices begin at
1163ab693e5aSLois Curfman McInnes    one, not zero.  See the users manual for further details.
1164ff756334SLois Curfman McInnes 
1165ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1166ff756334SLois Curfman McInnes    (possibly both).
1167ff756334SLois Curfman McInnes 
1168e0245417SLois Curfman McInnes    Storage Information:
1169e0245417SLois Curfman McInnes    For a square global matrix we define each processor's diagonal portion
1170e0245417SLois Curfman McInnes    to be its local rows and the corresponding columns (a square submatrix);
1171e0245417SLois Curfman McInnes    each processor's off-diagonal portion encompasses the remainder of the
1172e0245417SLois Curfman McInnes    local matrix (a rectangular submatrix).
1173e0245417SLois Curfman McInnes 
1174e0245417SLois Curfman McInnes    The user can specify preallocated storage for the diagonal part of
1175e0245417SLois Curfman McInnes    the local submatrix with either d_nz or d_nnz (not both). Set both
1176e0245417SLois Curfman McInnes    d_nz and d_nnz to zero for PETSc to control dynamic memory allocation.
1177e0245417SLois Curfman McInnes    Likewise, specify preallocated storage for the off-diagonal part of
1178e0245417SLois Curfman McInnes    the local submatrix with o_nz or o_nnz (not both).
1179ff756334SLois Curfman McInnes 
1180dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel
1181ff756334SLois Curfman McInnes 
1182dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatCreateSequentialAIJ(), MatSetValues()
11838a729477SBarry Smith @*/
11841eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,
11851eb62cbbSBarry Smith                  int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat)
11868a729477SBarry Smith {
11878a729477SBarry Smith   Mat          mat;
118844a69424SLois Curfman McInnes   Mat_MPIAIJ   *aij;
11896abc6512SBarry Smith   int          ierr, i,sum[2],work[2];
11908a729477SBarry Smith   *newmat         = 0;
119144a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm);
1192a5a9c739SBarry Smith   PLogObjectCreate(mat);
119378b31e54SBarry Smith   mat->data       = (void *) (aij = PETSCNEW(Mat_MPIAIJ)); CHKPTRQ(aij);
11948a729477SBarry Smith   mat->ops        = &MatOps;
119544a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
119644a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
11978a729477SBarry Smith   mat->factor     = 0;
1198d6dfbf8fSBarry Smith 
119907a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
12001eb62cbbSBarry Smith   MPI_Comm_rank(comm,&aij->mytid);
12011eb62cbbSBarry Smith   MPI_Comm_size(comm,&aij->numtids);
12021eb62cbbSBarry Smith 
1203dbd7a890SLois Curfman McInnes   if (M == PETSC_DECIDE || N == PETSC_DECIDE) {
12041eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1205d6dfbf8fSBarry Smith     MPI_Allreduce((void *) work,(void *) sum,2,MPI_INT,MPI_SUM,comm );
1206dbd7a890SLois Curfman McInnes     if (M == PETSC_DECIDE) M = sum[0];
1207dbd7a890SLois Curfman McInnes     if (N == PETSC_DECIDE) N = sum[1];
12081eb62cbbSBarry Smith   }
1209dbd7a890SLois Curfman McInnes   if (m == PETSC_DECIDE)
1210dbd7a890SLois Curfman McInnes     {m = M/aij->numtids + ((M % aij->numtids) > aij->mytid);}
1211dbd7a890SLois Curfman McInnes   if (n == PETSC_DECIDE)
1212dbd7a890SLois Curfman McInnes     {n = N/aij->numtids + ((N % aij->numtids) > aij->mytid);}
12138a729477SBarry Smith   aij->m = m;
12148a729477SBarry Smith   aij->n = n;
12151eb62cbbSBarry Smith   aij->N = N;
12161eb62cbbSBarry Smith   aij->M = M;
12171eb62cbbSBarry Smith 
12181eb62cbbSBarry Smith   /* build local table of row and column ownerships */
121978b31e54SBarry Smith   aij->rowners = (int *) PETSCMALLOC(2*(aij->numtids+2)*sizeof(int));
122078b31e54SBarry Smith   CHKPTRQ(aij->rowners);
1221464493b3SBarry Smith   PLogObjectMemory(mat,2*(aij->numtids+2)*sizeof(int)+sizeof(struct _Mat)+
1222464493b3SBarry Smith                        sizeof(Mat_MPIAIJ));
12231eb62cbbSBarry Smith   aij->cowners = aij->rowners + aij->numtids + 1;
12241eb62cbbSBarry Smith   MPI_Allgather(&m,1,MPI_INT,aij->rowners+1,1,MPI_INT,comm);
12251eb62cbbSBarry Smith   aij->rowners[0] = 0;
12261eb62cbbSBarry Smith   for ( i=2; i<=aij->numtids; i++ ) {
12271eb62cbbSBarry Smith     aij->rowners[i] += aij->rowners[i-1];
12288a729477SBarry Smith   }
12291eb62cbbSBarry Smith   aij->rstart = aij->rowners[aij->mytid];
12301eb62cbbSBarry Smith   aij->rend   = aij->rowners[aij->mytid+1];
12311eb62cbbSBarry Smith   MPI_Allgather(&n,1,MPI_INT,aij->cowners+1,1,MPI_INT,comm);
12321eb62cbbSBarry Smith   aij->cowners[0] = 0;
12331eb62cbbSBarry Smith   for ( i=2; i<=aij->numtids; i++ ) {
12341eb62cbbSBarry Smith     aij->cowners[i] += aij->cowners[i-1];
12358a729477SBarry Smith   }
12361eb62cbbSBarry Smith   aij->cstart = aij->cowners[aij->mytid];
12371eb62cbbSBarry Smith   aij->cend   = aij->cowners[aij->mytid+1];
12388a729477SBarry Smith 
12398a729477SBarry Smith 
12406b5873e3SBarry Smith   ierr = MatCreateSequentialAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&aij->A);
124178b31e54SBarry Smith   CHKERRQ(ierr);
1242a5a9c739SBarry Smith   PLogObjectParent(mat,aij->A);
12436b5873e3SBarry Smith   ierr = MatCreateSequentialAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&aij->B);
124478b31e54SBarry Smith   CHKERRQ(ierr);
1245a5a9c739SBarry Smith   PLogObjectParent(mat,aij->B);
12468a729477SBarry Smith 
12471eb62cbbSBarry Smith   /* build cache for off array entries formed */
124878b31e54SBarry Smith   ierr = StashBuild_Private(&aij->stash); CHKERRQ(ierr);
12499e25ed09SBarry Smith   aij->colmap    = 0;
12509e25ed09SBarry Smith   aij->garray    = 0;
12518a729477SBarry Smith 
12521eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
12531eb62cbbSBarry Smith   aij->lvec      = 0;
12541eb62cbbSBarry Smith   aij->Mvctx     = 0;
1255d6dfbf8fSBarry Smith   aij->assembled = 0;
12568a729477SBarry Smith 
1257d6dfbf8fSBarry Smith   *newmat = mat;
1258d6dfbf8fSBarry Smith   return 0;
1259d6dfbf8fSBarry Smith }
1260c74985f6SBarry Smith 
1261ff7509f2SBarry Smith static int MatCopyPrivate_MPIAIJ(Mat matin,Mat *newmat)
1262d6dfbf8fSBarry Smith {
1263d6dfbf8fSBarry Smith   Mat        mat;
126444a69424SLois Curfman McInnes   Mat_MPIAIJ *aij,*oldmat = (Mat_MPIAIJ *) matin->data;
12652ee70a88SLois Curfman McInnes   int        ierr, len;
1266d6dfbf8fSBarry Smith   *newmat      = 0;
1267d6dfbf8fSBarry Smith 
1268bbb6d6a8SBarry Smith   if (!oldmat->assembled)
1269bbb6d6a8SBarry Smith     SETERRQ(1,"MatCopyPrivate_MPIAIJ:Cannot copy unassembled matrix");
127044a69424SLois Curfman McInnes   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm);
1271a5a9c739SBarry Smith   PLogObjectCreate(mat);
127278b31e54SBarry Smith   mat->data       = (void *) (aij = PETSCNEW(Mat_MPIAIJ)); CHKPTRQ(aij);
1273d6dfbf8fSBarry Smith   mat->ops        = &MatOps;
127444a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
127544a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
1276d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1277d6dfbf8fSBarry Smith 
1278d6dfbf8fSBarry Smith   aij->m          = oldmat->m;
1279d6dfbf8fSBarry Smith   aij->n          = oldmat->n;
1280d6dfbf8fSBarry Smith   aij->M          = oldmat->M;
1281d6dfbf8fSBarry Smith   aij->N          = oldmat->N;
1282d6dfbf8fSBarry Smith 
1283d6dfbf8fSBarry Smith   aij->assembled  = 1;
1284d6dfbf8fSBarry Smith   aij->rstart     = oldmat->rstart;
1285d6dfbf8fSBarry Smith   aij->rend       = oldmat->rend;
1286d6dfbf8fSBarry Smith   aij->cstart     = oldmat->cstart;
1287d6dfbf8fSBarry Smith   aij->cend       = oldmat->cend;
1288d6dfbf8fSBarry Smith   aij->numtids    = oldmat->numtids;
1289d6dfbf8fSBarry Smith   aij->mytid      = oldmat->mytid;
129007a0e7adSLois Curfman McInnes   aij->insertmode = NOTSETVALUES;
1291d6dfbf8fSBarry Smith 
129278b31e54SBarry Smith   aij->rowners        = (int *) PETSCMALLOC( (aij->numtids+1)*sizeof(int) );
129378b31e54SBarry Smith   CHKPTRQ(aij->rowners);
1294464493b3SBarry Smith   PLogObjectMemory(mat,(aij->numtids+1)*sizeof(int)+sizeof(struct _Mat)+
1295464493b3SBarry Smith                        sizeof(Mat_MPIAIJ));
129678b31e54SBarry Smith   PETSCMEMCPY(aij->rowners,oldmat->rowners,(aij->numtids+1)*sizeof(int));
129778b31e54SBarry Smith   ierr = StashInitialize_Private(&aij->stash); CHKERRQ(ierr);
12982ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
129978b31e54SBarry Smith     aij->colmap      = (int *) PETSCMALLOC( (aij->N)*sizeof(int) );
130078b31e54SBarry Smith     CHKPTRQ(aij->colmap);
1301464493b3SBarry Smith     PLogObjectMemory(mat,(aij->N)*sizeof(int));
130278b31e54SBarry Smith     PETSCMEMCPY(aij->colmap,oldmat->colmap,(aij->N)*sizeof(int));
13032ee70a88SLois Curfman McInnes   } else aij->colmap = 0;
13042ee70a88SLois Curfman McInnes   if (oldmat->garray && (len = ((Mat_AIJ *) (oldmat->B->data))->n)) {
130578b31e54SBarry Smith     aij->garray      = (int *) PETSCMALLOC(len*sizeof(int) ); CHKPTRQ(aij->garray);
1306464493b3SBarry Smith     PLogObjectMemory(mat,len*sizeof(int));
130778b31e54SBarry Smith     PETSCMEMCPY(aij->garray,oldmat->garray,len*sizeof(int));
13082ee70a88SLois Curfman McInnes   } else aij->garray = 0;
1309d6dfbf8fSBarry Smith 
131078b31e54SBarry Smith   ierr =  VecDuplicate(oldmat->lvec,&aij->lvec); CHKERRQ(ierr);
1311a5a9c739SBarry Smith   PLogObjectParent(mat,aij->lvec);
131278b31e54SBarry Smith   ierr =  VecScatterCtxCopy(oldmat->Mvctx,&aij->Mvctx); CHKERRQ(ierr);
1313a5a9c739SBarry Smith   PLogObjectParent(mat,aij->Mvctx);
131478b31e54SBarry Smith   ierr =  MatConvert(oldmat->A,MATSAME,&aij->A); CHKERRQ(ierr);
1315a5a9c739SBarry Smith   PLogObjectParent(mat,aij->A);
131678b31e54SBarry Smith   ierr =  MatConvert(oldmat->B,MATSAME,&aij->B); CHKERRQ(ierr);
1317a5a9c739SBarry Smith   PLogObjectParent(mat,aij->B);
13188a729477SBarry Smith   *newmat = mat;
13198a729477SBarry Smith   return 0;
13208a729477SBarry Smith }
1321