xref: /petsc/src/ksp/pc/impls/redundant/redundant.c (revision ca6161988a13f5f1c12dfd3b7a0b9d56a75fe920)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
43f457be1SHong Zhang   This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
54b9ad928SBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
74b9ad928SBarry Smith #include "petscksp.h"
84b9ad928SBarry Smith 
94b9ad928SBarry Smith typedef struct {
104b9ad928SBarry Smith   PC         pc;                   /* actual preconditioner used on each processor */
113f457be1SHong Zhang   Vec        xsub,ysub;            /* vectors of a subcommunicator to hold parallel vectors of pc->comm */
123f457be1SHong Zhang   Vec        xdup,ydup;            /* parallel vector that congregates xsub or ysub facilitating vector scattering */
134b9ad928SBarry Smith   Mat        *pmats;               /* matrix and optional preconditioner matrix */
143f457be1SHong Zhang   Mat        pmats_sub;            /* matrix and optional preconditioner matrix */
153f457be1SHong Zhang   VecScatter scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */
164b9ad928SBarry Smith   PetscTruth useparallelmat;
173f457be1SHong Zhang   MPI_Comm   subcomm;              /* processors belong to a subcommunicator implement a PC in parallel */
183f457be1SHong Zhang   MPI_Comm   dupcomm;              /* processors belong to pc->comm with their rank remapped in the way
193f457be1SHong Zhang                                       that vector xdup/ydup has contiguous rank while appending xsub/ysub along their colors */
203f457be1SHong Zhang   PetscInt   nsubcomm;             /* num of subcommunicators, which equals the num of redundant matrix systems */
213f457be1SHong Zhang   PetscInt   color;                /* color of processors in a subcommunicator */
224b9ad928SBarry Smith } PC_Redundant;
234b9ad928SBarry Smith 
244b9ad928SBarry Smith #undef __FUNCT__
254b9ad928SBarry Smith #define __FUNCT__ "PCView_Redundant"
266849ba73SBarry Smith static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer)
274b9ad928SBarry Smith {
284b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
29dfbe8321SBarry Smith   PetscErrorCode ierr;
3013f74950SBarry Smith   PetscMPIInt    rank;
3132077d6dSBarry Smith   PetscTruth     iascii,isstring;
324b9ad928SBarry Smith   PetscViewer    sviewer;
334b9ad928SBarry Smith 
344b9ad928SBarry Smith   PetscFunctionBegin;
354b9ad928SBarry Smith   ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr);
3632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
374b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3832077d6dSBarry Smith   if (iascii) {
394b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Redundant solver preconditioner: Actual PC follows\n");CHKERRQ(ierr);
404b9ad928SBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
414b9ad928SBarry Smith     if (!rank) {
424b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
434b9ad928SBarry Smith       ierr = PCView(red->pc,sviewer);CHKERRQ(ierr);
444b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
454b9ad928SBarry Smith     }
464b9ad928SBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
474b9ad928SBarry Smith   } else if (isstring) {
484b9ad928SBarry Smith     ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr);
494b9ad928SBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
504b9ad928SBarry Smith     if (!rank) {
514b9ad928SBarry Smith       ierr = PCView(red->pc,sviewer);CHKERRQ(ierr);
524b9ad928SBarry Smith     }
534b9ad928SBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
544b9ad928SBarry Smith   } else {
5579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name);
564b9ad928SBarry Smith   }
574b9ad928SBarry Smith   PetscFunctionReturn(0);
584b9ad928SBarry Smith }
594b9ad928SBarry Smith 
603f457be1SHong Zhang #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
613f457be1SHong Zhang #include "private/vecimpl.h"
623f457be1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"   /*I "petscmat.h" I*/
633f457be1SHong Zhang #include "src/mat/impls/aij/seq/aij.h"      /*I "petscmat.h" I*/
643f457be1SHong Zhang 
653f457be1SHong Zhang #undef __FUNCT__
663f457be1SHong Zhang #define __FUNCT__ "MatGetRedundantMatrix"
673f457be1SHong Zhang PetscErrorCode MatGetRedundantMatrix_AIJ(PC pc,Mat mat,MPI_Comm subcomm,PetscInt mlocal_sub,Mat *matredundant)
683f457be1SHong Zhang {
693f457be1SHong Zhang   PetscMPIInt    rank,size,subrank,subsize;
703f457be1SHong Zhang   MPI_Comm       comm=mat->comm;
713f457be1SHong Zhang   PetscErrorCode ierr;
723f457be1SHong Zhang   PC_Redundant   *red=(PC_Redundant*)pc->data;
733f457be1SHong Zhang   PetscInt       nsubcomm=red->nsubcomm,nsends,nrecvs,i,prid=100,itmp;
743f457be1SHong Zhang   PetscMPIInt    *send_rank,*recv_rank;
75*ca616198SHong Zhang   PetscInt       *rowrange=pc->pmat->rmap.range,nzlocal;
763f457be1SHong Zhang   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)mat->data;
773f457be1SHong Zhang   Mat            A=aij->A,B=aij->B;
783f457be1SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data;
793f457be1SHong Zhang   Mat            C;
803f457be1SHong Zhang 
813f457be1SHong Zhang   PetscInt       nleftover,np_subcomm,j;
823f457be1SHong Zhang   PetscInt       nz_A,nz_B,*sbuf_j;
833f457be1SHong Zhang   PetscScalar    *sbuf_a;
843f457be1SHong Zhang   PetscInt       cstart=mat->cmap.rstart,cend=mat->cmap.rend,row,nzA,nzB,ncols,*cworkA,*cworkB;
853f457be1SHong Zhang   PetscInt       rstart=mat->rmap.rstart,rend=mat->rmap.rend,*bmap=aij->garray,M,N;
863f457be1SHong Zhang   PetscInt       *cols,ctmp,lwrite,*rptr,l;
873f457be1SHong Zhang   PetscScalar    *vals,*aworkA,*aworkB;
883f457be1SHong Zhang 
893f457be1SHong Zhang   PetscMPIInt    tag1,tag2,tag3,imdex;
903f457be1SHong Zhang   MPI_Request    *s_waits1,*s_waits2,*s_waits3,*r_waits1,*r_waits2,*r_waits3;
913f457be1SHong Zhang   MPI_Status     recv_status,*send_status;
923f457be1SHong Zhang   PetscInt       *sbuf_nz,*rbuf_nz,count;
933f457be1SHong Zhang 
943f457be1SHong Zhang   PetscInt     **rbuf_j;
953f457be1SHong Zhang   PetscScalar  **rbuf_a;
963f457be1SHong Zhang 
973f457be1SHong Zhang   PetscFunctionBegin;
983f457be1SHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-prid",&prid,PETSC_NULL);CHKERRQ(ierr);
993f457be1SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1003f457be1SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1013f457be1SHong Zhang   ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr);
1023f457be1SHong Zhang   ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
1033f457be1SHong Zhang   /*
1043f457be1SHong Zhang   ierr = PetscSynchronizedPrintf(comm, "[%d] subrank %d, subsize %d\n",rank,subrank,subsize);
1053f457be1SHong Zhang   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
1063f457be1SHong Zhang   */
1073f457be1SHong Zhang   /* get the destination processors */
1083f457be1SHong Zhang   ierr = PetscMalloc((2*size+1)*sizeof(PetscMPIInt),&send_rank);
1093f457be1SHong Zhang   recv_rank = send_rank + size;
1103f457be1SHong Zhang   np_subcomm = size/nsubcomm;
1113f457be1SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
1123f457be1SHong Zhang   nsends = 0; nrecvs = 0;
1133f457be1SHong Zhang   for (i=0; i<size; i++){ /* i=rank*/
1143f457be1SHong Zhang     if (subrank == i/nsubcomm && rank != i){ /* my_subrank == other's subrank */
1153f457be1SHong Zhang       send_rank[nsends] = i; nsends++;
1163f457be1SHong Zhang       recv_rank[nrecvs++] = i;
1173f457be1SHong Zhang     }
1183f457be1SHong Zhang   }
1193f457be1SHong Zhang   if (rank >= size - nleftover){/* this proc is a leftover processor */
1203f457be1SHong Zhang     i = size-nleftover-1;
1213f457be1SHong Zhang     j = 0;
1223f457be1SHong Zhang     while (j < nsubcomm - nleftover){
1233f457be1SHong Zhang       send_rank[nsends++] = i;
1243f457be1SHong Zhang       i--; j++;
1253f457be1SHong Zhang     }
1263f457be1SHong Zhang   }
1273f457be1SHong Zhang 
1283f457be1SHong Zhang   if (nleftover && subsize == size/nsubcomm && subrank==subsize-1){ /* this proc recvs from leftover processors */
1293f457be1SHong Zhang     for (i=0; i<nleftover; i++){
1303f457be1SHong Zhang       recv_rank[nrecvs++] = size-nleftover+i;
1313f457be1SHong Zhang     }
1323f457be1SHong Zhang   }
1333f457be1SHong Zhang   if (prid == rank){
1343f457be1SHong Zhang     printf("[%d] sends to ",rank);
1353f457be1SHong Zhang     for (i=0; i<nsends; i++) printf(" [%d],",send_rank[i]);
1363f457be1SHong Zhang     printf("  \n");
1373f457be1SHong Zhang   }
1383f457be1SHong Zhang   /*
1393f457be1SHong Zhang   ierr = PetscSynchronizedPrintf(comm, "[%d] nsends %d, nrecvs %d\n",rank,nsends,nrecvs);
1403f457be1SHong Zhang   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
1413f457be1SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1423f457be1SHong Zhang   */
1433f457be1SHong Zhang 
1443f457be1SHong Zhang   /* get this processor's nzlocal=nz_A+nz_B */
1453f457be1SHong Zhang   nz_A = a->nz; nz_B = b->nz;
1463f457be1SHong Zhang   nzlocal = nz_A + nz_B;
1473f457be1SHong Zhang 
1483f457be1SHong Zhang   /* allocate sbuf_j, sbuf_a, then copy mat's local entries into the buffers */
1493f457be1SHong Zhang   itmp = nzlocal + rowrange[rank+1] - rowrange[rank] + 2;
1503f457be1SHong Zhang   ierr = PetscMalloc(itmp*sizeof(PetscInt),&sbuf_j);CHKERRQ(ierr);
1513f457be1SHong Zhang   ierr = PetscMalloc((nzlocal+1)*sizeof(PetscScalar),&sbuf_a);CHKERRQ(ierr);
1523f457be1SHong Zhang 
1533f457be1SHong Zhang     rptr = sbuf_j;
1543f457be1SHong Zhang     cols = sbuf_j + rend-rstart + 1;
1553f457be1SHong Zhang     vals = sbuf_a;
1563f457be1SHong Zhang     rptr[0] = 0;
1573f457be1SHong Zhang     for (i=0; i<rend-rstart; i++){
1583f457be1SHong Zhang       row = i + rstart;
1593f457be1SHong Zhang       nzA    = a->i[i+1] - a->i[i]; nzB = b->i[i+1] - b->i[i];
1603f457be1SHong Zhang       ncols  = nzA + nzB;
1613f457be1SHong Zhang       cworkA = a->j + a->i[i]; cworkB = b->j + b->i[i];
1623f457be1SHong Zhang       aworkA = a->a + a->i[i]; aworkB = b->a + b->i[i];
1633f457be1SHong Zhang       /* load the column indices for this row into cols */
1643f457be1SHong Zhang       lwrite = 0;
1653f457be1SHong Zhang       for (l=0; l<nzB; l++) {
1663f457be1SHong Zhang         if ((ctmp = bmap[cworkB[l]]) < cstart){
1673f457be1SHong Zhang           vals[lwrite]   = aworkB[l];
1683f457be1SHong Zhang           cols[lwrite++] = ctmp;
1693f457be1SHong Zhang         }
1703f457be1SHong Zhang       }
1713f457be1SHong Zhang       for (l=0; l<nzA; l++){
1723f457be1SHong Zhang         vals[lwrite]   = aworkA[l];
1733f457be1SHong Zhang         cols[lwrite++] = cstart + cworkA[l];
1743f457be1SHong Zhang       }
1753f457be1SHong Zhang       for (l=0; l<nzB; l++) {
1763f457be1SHong Zhang         if ((ctmp = bmap[cworkB[l]]) >= cend){
1773f457be1SHong Zhang           vals[lwrite]   = aworkB[l];
1783f457be1SHong Zhang           cols[lwrite++] = ctmp;
1793f457be1SHong Zhang         }
1803f457be1SHong Zhang       }
1813f457be1SHong Zhang       /* insert local matrix values into C */
182*ca616198SHong Zhang       /* ierr = MatSetValues(C,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr); */
1833f457be1SHong Zhang 
1843f457be1SHong Zhang       vals += ncols;
1853f457be1SHong Zhang       cols += ncols;
1863f457be1SHong Zhang       rptr[i+1] = rptr[i] + ncols;
1873f457be1SHong Zhang     }
1883f457be1SHong Zhang     /*
1893f457be1SHong Zhang     ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1903f457be1SHong Zhang     ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1913f457be1SHong Zhang     */
1923f457be1SHong Zhang     if (rptr[rend-rstart] != a->nz + b->nz) SETERRQ4(1, "rptr[%d] %d != %d + %d",rend-rstart,rptr[rend-rstart+1],a->nz,b->nz);
1933f457be1SHong Zhang 
1943f457be1SHong Zhang   /* send nzlocal to others, and recv other's nzlocal */
1953f457be1SHong Zhang   /*--------------------------------------------------*/
1963f457be1SHong Zhang   ierr = PetscMalloc3(nsends+nrecvs+1,PetscInt,&sbuf_nz,nrecvs,PetscInt*,&rbuf_j,nrecvs,PetscScalar*,&rbuf_a);CHKERRQ(ierr);
1973f457be1SHong Zhang   rbuf_nz = sbuf_nz + nsends;
1983f457be1SHong Zhang 
1993f457be1SHong Zhang   ierr = PetscMalloc2(3*(nsends + nrecvs)+1,MPI_Request,&s_waits1,nsends+1,MPI_Status,&send_status);CHKERRQ(ierr);
2003f457be1SHong Zhang   s_waits2 = s_waits1 + nsends;
2013f457be1SHong Zhang   s_waits3 = s_waits2 + nsends;
2023f457be1SHong Zhang   r_waits1 = s_waits3 + nsends;
2033f457be1SHong Zhang   r_waits2 = r_waits1 + nrecvs;
2043f457be1SHong Zhang   r_waits3 = r_waits2 + nrecvs;
2053f457be1SHong Zhang 
2063f457be1SHong Zhang   /* get some new tags to keep the communication clean */
2073f457be1SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)A,&tag1);CHKERRQ(ierr);
2083f457be1SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)A,&tag2);CHKERRQ(ierr);
2093f457be1SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)A,&tag3);CHKERRQ(ierr);
2103f457be1SHong Zhang 
2113f457be1SHong Zhang   /* post receives of other's nzlocal */
2123f457be1SHong Zhang   for (i=0; i<nrecvs; i++){
2133f457be1SHong Zhang     ierr = MPI_Irecv(rbuf_nz+i,1,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,r_waits1+i);CHKERRQ(ierr);
2143f457be1SHong Zhang   }
2153f457be1SHong Zhang 
2163f457be1SHong Zhang   /* send nzlocal to others */
2173f457be1SHong Zhang   for (i=0; i<nsends; i++){
2183f457be1SHong Zhang     sbuf_nz[i] = nzlocal;
2193f457be1SHong Zhang     ierr = MPI_Isend(sbuf_nz+i,1,MPIU_INT,send_rank[i],tag1,comm,s_waits1+i);CHKERRQ(ierr);
2203f457be1SHong Zhang     if (prid == rank) printf(" [%d] sends nz %d to [%d]\n",rank,nzlocal,send_rank[i]);
2213f457be1SHong Zhang   }
2223f457be1SHong Zhang 
2233f457be1SHong Zhang   /* wait on receives of nzlocal; allocate space for rbuf_j, rbuf_a */
2243f457be1SHong Zhang   count = nrecvs;
2253f457be1SHong Zhang   while (count) {
2263f457be1SHong Zhang     ierr = MPI_Waitany(nrecvs,r_waits1,&imdex,&recv_status);CHKERRQ(ierr);
2273f457be1SHong Zhang     recv_rank[imdex] = recv_status.MPI_SOURCE;
2283f457be1SHong Zhang     /* allocate rbuf_a and rbuf_j; then post receives of rbuf_a and rbuf_j */
2293f457be1SHong Zhang     ierr = PetscMalloc((rbuf_nz[imdex]+1)*sizeof(PetscScalar),&rbuf_a[imdex]);CHKERRQ(ierr);
2303f457be1SHong Zhang     ierr = MPI_Irecv(rbuf_a[imdex],rbuf_nz[imdex],MPIU_SCALAR,recv_status.MPI_SOURCE,tag3,comm,r_waits3+imdex);CHKERRQ(ierr);
2313f457be1SHong Zhang 
2323f457be1SHong Zhang     itmp = rowrange[recv_status.MPI_SOURCE+1] - rowrange[recv_status.MPI_SOURCE]; /* number of expected mat->i */
2333f457be1SHong Zhang     rbuf_nz[imdex] += itmp+2;
2343f457be1SHong Zhang     ierr = PetscMalloc(rbuf_nz[imdex]*sizeof(PetscInt),&rbuf_j[imdex]);CHKERRQ(ierr);
2353f457be1SHong Zhang     ierr = MPI_Irecv(rbuf_j[imdex],rbuf_nz[imdex],MPIU_INT,recv_status.MPI_SOURCE,tag2,comm,r_waits2+imdex);CHKERRQ(ierr);
2363f457be1SHong Zhang     count--;
2373f457be1SHong Zhang   }
2383f457be1SHong Zhang 
2393f457be1SHong Zhang   /* wait on sends of nzlocal */
2403f457be1SHong Zhang   if (nsends) {ierr = MPI_Waitall(nsends,s_waits1,send_status);CHKERRQ(ierr);}
2413f457be1SHong Zhang 
2423f457be1SHong Zhang   /* send mat->i,j and mat->a to others, and recv from other's */
2433f457be1SHong Zhang   /*-----------------------------------------------------------*/
2443f457be1SHong Zhang   for (i=0; i<nsends; i++){
2453f457be1SHong Zhang     ierr = MPI_Isend(sbuf_a,nzlocal,MPIU_SCALAR,send_rank[i],tag3,comm,s_waits3+i);CHKERRQ(ierr);
2463f457be1SHong Zhang     itmp = nzlocal + rowrange[rank+1] - rowrange[rank] + 1;
2473f457be1SHong Zhang     ierr = MPI_Isend(sbuf_j,itmp,MPIU_INT,send_rank[i],tag2,comm,s_waits2+i);CHKERRQ(ierr);
2483f457be1SHong Zhang   }
2493f457be1SHong Zhang 
2503f457be1SHong Zhang   /* wait on receives of mat->i,j and mat->a */
2513f457be1SHong Zhang   /*-----------------------------------------*/
2523f457be1SHong Zhang   count = nrecvs;
2533f457be1SHong Zhang   while (count) {
2543f457be1SHong Zhang     ierr = MPI_Waitany(nrecvs,r_waits3,&imdex,&recv_status);CHKERRQ(ierr);
2553f457be1SHong Zhang     if (recv_rank[imdex] != recv_status.MPI_SOURCE) SETERRQ2(1, "recv_rank %d != MPI_SOURCE %d",recv_rank[imdex],recv_status.MPI_SOURCE);
2563f457be1SHong Zhang     ierr = MPI_Waitany(nrecvs,r_waits2,&imdex,&recv_status);CHKERRQ(ierr);
2573f457be1SHong Zhang     count--;
2583f457be1SHong Zhang   }
2593f457be1SHong Zhang 
2603f457be1SHong Zhang   /* wait on sends of mat->i,j and mat->a */
2613f457be1SHong Zhang   /*--------------------------------------*/
2623f457be1SHong Zhang   if (nsends) {
2633f457be1SHong Zhang     ierr = MPI_Waitall(nsends,s_waits3,send_status);CHKERRQ(ierr);
2643f457be1SHong Zhang     ierr = MPI_Waitall(nsends,s_waits2,send_status);CHKERRQ(ierr);
2653f457be1SHong Zhang   }
2663f457be1SHong Zhang   ierr = PetscFree(sbuf_nz);CHKERRQ(ierr);
2673f457be1SHong Zhang   ierr = PetscFree2(s_waits1,send_status);CHKERRQ(ierr);
2683f457be1SHong Zhang 
2693f457be1SHong Zhang   /* create redundant matrix */
2703f457be1SHong Zhang   /*-------------------------*/
2713f457be1SHong Zhang   ierr = MatCreate(subcomm,&C);CHKERRQ(ierr);
2723f457be1SHong Zhang   ierr = MatSetSizes(C,mlocal_sub,mlocal_sub,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
2733f457be1SHong Zhang   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
2743f457be1SHong Zhang 
2753f457be1SHong Zhang   /* insert local matrix entries */
2763f457be1SHong Zhang   rptr = sbuf_j;
2773f457be1SHong Zhang   cols = sbuf_j + rend-rstart + 1;
2783f457be1SHong Zhang   vals = sbuf_a;
2793f457be1SHong Zhang   for (i=0; i<rend-rstart; i++){
2803f457be1SHong Zhang     row   = i + rstart;
2813f457be1SHong Zhang     ncols = rptr[i+1] - rptr[i];
2823f457be1SHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
2833f457be1SHong Zhang     vals += ncols;
2843f457be1SHong Zhang     cols += ncols;
2853f457be1SHong Zhang   }
2863f457be1SHong Zhang   /* insert received matrix entries */
2873f457be1SHong Zhang   for (imdex=0; imdex<nrecvs; imdex++){
2883f457be1SHong Zhang 
2893f457be1SHong Zhang     rstart = rowrange[recv_rank[imdex]];
2903f457be1SHong Zhang     rend   = rowrange[recv_rank[imdex]+1];
2913f457be1SHong Zhang     rptr = rbuf_j[imdex];
2923f457be1SHong Zhang     cols = rbuf_j[imdex] + rend-rstart + 1;
2933f457be1SHong Zhang     vals = rbuf_a[imdex];
2943f457be1SHong Zhang     for (i=0; i<rend-rstart; i++){
2953f457be1SHong Zhang       row   = i + rstart;
2963f457be1SHong Zhang       ncols = rptr[i+1] - rptr[i];
2973f457be1SHong Zhang       ierr = MatSetValues(C,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
2983f457be1SHong Zhang       vals += ncols;
2993f457be1SHong Zhang       cols += ncols;
3003f457be1SHong Zhang     }
3013f457be1SHong Zhang   }
3023f457be1SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3033f457be1SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3043f457be1SHong Zhang   ierr = MatGetSize(C,&M,&N);CHKERRQ(ierr);
3053f457be1SHong Zhang   if (M != mat->rmap.N || N != mat->cmap.N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"redundant mat size %d != input mat size %d",M,mat->rmap.N);
3063f457be1SHong Zhang   *matredundant = C;
3073f457be1SHong Zhang 
3083f457be1SHong Zhang   /* free space */
3093f457be1SHong Zhang   ierr = PetscFree(send_rank);CHKERRQ(ierr);
3103f457be1SHong Zhang   ierr = PetscFree(sbuf_j);CHKERRQ(ierr);
3113f457be1SHong Zhang   ierr = PetscFree(sbuf_a);CHKERRQ(ierr);
3123f457be1SHong Zhang   for (i=0; i<nrecvs; i++){
3133f457be1SHong Zhang     ierr = PetscFree(rbuf_j[i]);CHKERRQ(ierr);
3143f457be1SHong Zhang     ierr = PetscFree(rbuf_a[i]);CHKERRQ(ierr);
3153f457be1SHong Zhang   }
3163f457be1SHong Zhang   ierr = PetscFree3(sbuf_nz,rbuf_j,rbuf_a);CHKERRQ(ierr);
3173f457be1SHong Zhang   PetscFunctionReturn(0);
3183f457be1SHong Zhang }
3193f457be1SHong Zhang 
3204b9ad928SBarry Smith #undef __FUNCT__
3214b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant"
3226849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc)
3234b9ad928SBarry Smith {
3244b9ad928SBarry Smith   PC_Redundant   *red  = (PC_Redundant*)pc->data;
325dfbe8321SBarry Smith   PetscErrorCode ierr;
3263f457be1SHong Zhang   PetscInt       mstart,mend,mlocal,m;
32713f74950SBarry Smith   PetscMPIInt    size;
3284b9ad928SBarry Smith   IS             isl;
3294b9ad928SBarry Smith   MatReuse       reuse = MAT_INITIAL_MATRIX;
3304b9ad928SBarry Smith   MatStructure   str   = DIFFERENT_NONZERO_PATTERN;
3314b9ad928SBarry Smith   MPI_Comm       comm;
33223ce1328SBarry Smith   Vec            vec;
3334b9ad928SBarry Smith 
3343f457be1SHong Zhang   PetscInt    mlocal_sub;
3353f457be1SHong Zhang   PetscMPIInt subsize,subrank;
3363f457be1SHong Zhang   PetscInt    rstart_sub,rend_sub,mloc_sub;
3373f457be1SHong Zhang   Mat         Aredundant;
3383f457be1SHong Zhang 
3394b9ad928SBarry Smith   PetscFunctionBegin;
3403f457be1SHong Zhang   ierr = PetscObjectGetComm((PetscObject)pc->pmat,&comm);CHKERRQ(ierr);
34123ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
3424b9ad928SBarry Smith   ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr);
34323ce1328SBarry Smith   ierr = VecGetSize(vec,&m);CHKERRQ(ierr);
3444b9ad928SBarry Smith   if (!pc->setupcalled) {
3453f457be1SHong Zhang     /* create working vectors xsub/ysub and xdup/ydup */
34623ce1328SBarry Smith     ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr);
3473f457be1SHong Zhang     ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr);
3484b9ad928SBarry Smith 
3493f457be1SHong Zhang     /* get local size of xsub/ysub */
3503f457be1SHong Zhang     ierr = MPI_Comm_size(red->subcomm,&subsize);CHKERRQ(ierr);
3513f457be1SHong Zhang     ierr = MPI_Comm_rank(red->subcomm,&subrank);CHKERRQ(ierr);
3523f457be1SHong Zhang     rstart_sub = pc->pmat->rmap.range[red->nsubcomm*subrank]; /* rstart in xsub/ysub */
3533f457be1SHong Zhang     if (subrank+1 < subsize){
3543f457be1SHong Zhang       rend_sub = pc->pmat->rmap.range[red->nsubcomm*(subrank+1)];
3553f457be1SHong Zhang     } else {
3563f457be1SHong Zhang       rend_sub = m;
3573f457be1SHong Zhang     }
3583f457be1SHong Zhang     mloc_sub = rend_sub - rstart_sub;
3593f457be1SHong Zhang     ierr = VecCreateMPI(red->subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr);
3603f457be1SHong Zhang     /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
3613f457be1SHong Zhang     ierr = VecCreateMPIWithArray(red->subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr);
3623f457be1SHong Zhang 
3633f457be1SHong Zhang     /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it.
3643f457be1SHong Zhang        Note: we use communicator dupcomm, not pc->comm! */
3653f457be1SHong Zhang     ierr = VecCreateMPI(red->dupcomm,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr);
3663f457be1SHong Zhang     ierr = VecCreateMPIWithArray(red->dupcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr);
3673f457be1SHong Zhang 
3683f457be1SHong Zhang     /* create vec scatters */
3693f457be1SHong Zhang     if (!red->scatterin){
3703f457be1SHong Zhang       IS       is1,is2;
3713f457be1SHong Zhang       PetscInt *idx1,*idx2,i,j,k;
3723f457be1SHong Zhang       ierr = PetscMalloc(2*red->nsubcomm*mlocal*sizeof(PetscInt),&idx1);CHKERRQ(ierr);
3733f457be1SHong Zhang       idx2 = idx1 + red->nsubcomm*mlocal;
3743f457be1SHong Zhang       j = 0;
3753f457be1SHong Zhang       for (k=0; k<red->nsubcomm; k++){
3763f457be1SHong Zhang         for (i=mstart; i<mend; i++){
3773f457be1SHong Zhang           idx1[j]   = i;
3783f457be1SHong Zhang           idx2[j++] = i + m*k;
3793f457be1SHong Zhang         }
3803f457be1SHong Zhang       }
3813f457be1SHong Zhang       ierr = ISCreateGeneral(comm,red->nsubcomm*mlocal,idx1,&is1);CHKERRQ(ierr);
3823f457be1SHong Zhang       ierr = ISCreateGeneral(comm,red->nsubcomm*mlocal,idx2,&is2);CHKERRQ(ierr);
3833f457be1SHong Zhang       ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr);
3843f457be1SHong Zhang       ierr = ISDestroy(is1);CHKERRQ(ierr);
3853f457be1SHong Zhang       ierr = ISDestroy(is2);CHKERRQ(ierr);
3863f457be1SHong Zhang 
3873f457be1SHong Zhang       ierr = ISCreateStride(comm,mlocal,mstart+ red->color*m,1,&is1);CHKERRQ(ierr);
3883f457be1SHong Zhang       ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr);
3893f457be1SHong Zhang       ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr);
3903f457be1SHong Zhang       ierr = ISDestroy(is1);CHKERRQ(ierr);
3913f457be1SHong Zhang       ierr = ISDestroy(is2);CHKERRQ(ierr);
3923f457be1SHong Zhang       ierr = PetscFree(idx1);CHKERRQ(ierr);
3934b9ad928SBarry Smith     }
3944b9ad928SBarry Smith   }
39523ce1328SBarry Smith   ierr = VecDestroy(vec);CHKERRQ(ierr);
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
3983f457be1SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
3994b9ad928SBarry Smith   if (size == 1) {
4004b9ad928SBarry Smith     red->useparallelmat = PETSC_FALSE;
4014b9ad928SBarry Smith   }
4024b9ad928SBarry Smith 
4034b9ad928SBarry Smith   if (red->useparallelmat) {
4044b9ad928SBarry Smith     if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
4054b9ad928SBarry Smith       /* destroy old matrices */
4064b9ad928SBarry Smith       if (red->pmats) {
4074b9ad928SBarry Smith         ierr = MatDestroyMatrices(1,&red->pmats);CHKERRQ(ierr);
4084b9ad928SBarry Smith       }
4094b9ad928SBarry Smith     } else if (pc->setupcalled == 1) {
4104b9ad928SBarry Smith       reuse = MAT_REUSE_MATRIX;
4114b9ad928SBarry Smith       str   = SAME_NONZERO_PATTERN;
4124b9ad928SBarry Smith     }
4134b9ad928SBarry Smith 
4143f457be1SHong Zhang     /* ================== matrix ============= */
4153f457be1SHong Zhang     ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr);
4163f457be1SHong Zhang     ierr = MatGetRedundantMatrix_AIJ(pc,pc->pmat,red->subcomm,mlocal_sub,&Aredundant);CHKERRQ(ierr);
4173f457be1SHong Zhang 
4183f457be1SHong Zhang     /* grab the parallel matrix and put it into processors of a subcomminicator */
4194b9ad928SBarry Smith     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
4204b9ad928SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,1,&isl,&isl,reuse,&red->pmats);CHKERRQ(ierr);
4214b9ad928SBarry Smith     ierr = ISDestroy(isl);CHKERRQ(ierr);
4224b9ad928SBarry Smith 
4233f457be1SHong Zhang     /* ------- temporarily set a mpi matrix pmats_sub- provided by user! --*/
4243f457be1SHong Zhang     ierr = MatCreate(red->subcomm,&red->pmats_sub);CHKERRQ(ierr);
4253f457be1SHong Zhang     ierr = MatSetSizes(red->pmats_sub,mlocal_sub,mlocal_sub,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
4263f457be1SHong Zhang     ierr = MatSetFromOptions(red->pmats_sub);CHKERRQ(ierr);
4273f457be1SHong Zhang     {
4283f457be1SHong Zhang       PetscInt          Istart,Iend,ncols,i;
4293f457be1SHong Zhang       const PetscInt    *cols;
4303f457be1SHong Zhang       const PetscScalar *vals;
4313f457be1SHong Zhang       PetscTruth flg;
4323f457be1SHong Zhang       ierr = MatGetOwnershipRange(red->pmats_sub,&Istart,&Iend);CHKERRQ(ierr);
4333f457be1SHong Zhang       for (i=Istart; i<Iend; i++) {
4343f457be1SHong Zhang         ierr = MatGetRow(red->pmats[0],i,&ncols,&cols,&vals);CHKERRQ(ierr);
4353f457be1SHong Zhang         ierr = MatSetValues(red->pmats_sub,1,&i,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
4363f457be1SHong Zhang         ierr = MatRestoreRow(red->pmats[0],i,&ncols,&cols,&vals);CHKERRQ(ierr);
4373f457be1SHong Zhang       }
4383f457be1SHong Zhang       ierr = MatAssemblyBegin(red->pmats_sub,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4393f457be1SHong Zhang       ierr = MatAssemblyEnd(red->pmats_sub,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4403f457be1SHong Zhang 
4413f457be1SHong Zhang       ierr = MatEqual(red->pmats_sub,Aredundant,&flg);CHKERRQ(ierr);
4423f457be1SHong Zhang       if (!flg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"pmats_sub !=Aredundant ");
4433f457be1SHong Zhang       ierr = MatDestroy(red->pmats_sub);
4443f457be1SHong Zhang     }
4453f457be1SHong Zhang     red->pmats_sub = Aredundant;
4463f457be1SHong Zhang 
4473f457be1SHong Zhang     /* tell PC of the subcommunicator its operators */
4483f457be1SHong Zhang     ierr = PCSetOperators(red->pc,red->pmats_sub,red->pmats_sub,str);CHKERRQ(ierr);
4494b9ad928SBarry Smith   } else {
4504b9ad928SBarry Smith     ierr = PCSetOperators(red->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
4514b9ad928SBarry Smith   }
4524b9ad928SBarry Smith   ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr);
4534b9ad928SBarry Smith   ierr = PCSetUp(red->pc);CHKERRQ(ierr);
4544b9ad928SBarry Smith   PetscFunctionReturn(0);
4554b9ad928SBarry Smith }
4564b9ad928SBarry Smith 
4574b9ad928SBarry Smith #undef __FUNCT__
4584b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant"
4596849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
4604b9ad928SBarry Smith {
4614b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
462dfbe8321SBarry Smith   PetscErrorCode ierr;
4633f457be1SHong Zhang   PetscScalar    *array;
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith   PetscFunctionBegin;
4663f457be1SHong Zhang   /* scatter x to xdup */
4673f457be1SHong Zhang   ierr = VecScatterBegin(x,red->xdup,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr);
4683f457be1SHong Zhang   ierr = VecScatterEnd(x,red->xdup,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr);
4693f457be1SHong Zhang 
4703f457be1SHong Zhang   /* place xdup's local array into xsub */
4713f457be1SHong Zhang   ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr);
4723f457be1SHong Zhang   ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr);
4734b9ad928SBarry Smith 
4744b9ad928SBarry Smith   /* apply preconditioner on each processor */
4753f457be1SHong Zhang   ierr = PCApply(red->pc,red->xsub,red->ysub);CHKERRQ(ierr);
4763f457be1SHong Zhang   ierr = VecResetArray(red->xsub);CHKERRQ(ierr);
4773f457be1SHong Zhang   ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr);
4784b9ad928SBarry Smith 
4793f457be1SHong Zhang   /* place ysub's local array into ydup */
4803f457be1SHong Zhang   ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr);
4813f457be1SHong Zhang   ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr);
4823f457be1SHong Zhang 
4833f457be1SHong Zhang   /* scatter ydup to y */
4843f457be1SHong Zhang   ierr = VecScatterBegin(red->ydup,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr);
4853f457be1SHong Zhang   ierr = VecScatterEnd(red->ydup,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr);
4863f457be1SHong Zhang   ierr = VecResetArray(red->ydup);CHKERRQ(ierr);
4873f457be1SHong Zhang   ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr);
4884b9ad928SBarry Smith   PetscFunctionReturn(0);
4894b9ad928SBarry Smith }
4904b9ad928SBarry Smith 
4914b9ad928SBarry Smith #undef __FUNCT__
4924b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Redundant"
4936849ba73SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc)
4944b9ad928SBarry Smith {
4954b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
496dfbe8321SBarry Smith   PetscErrorCode ierr;
4974b9ad928SBarry Smith 
4984b9ad928SBarry Smith   PetscFunctionBegin;
4994b9ad928SBarry Smith   if (red->scatterin)  {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);}
5004b9ad928SBarry Smith   if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);}
5013f457be1SHong Zhang   if (red->ysub)       {ierr = VecDestroy(red->ysub);CHKERRQ(ierr);}
5023f457be1SHong Zhang   if (red->xsub)       {ierr = VecDestroy(red->xsub);CHKERRQ(ierr);}
5033f457be1SHong Zhang   if (red->xdup)       {ierr = VecDestroy(red->xdup);CHKERRQ(ierr);}
5043f457be1SHong Zhang   if (red->ydup)       {ierr = VecDestroy(red->ydup);CHKERRQ(ierr);}
5053f457be1SHong Zhang   if (red->pmats_sub) {
5063f457be1SHong Zhang     ierr = MatDestroy(red->pmats_sub);CHKERRQ(ierr);
5073f457be1SHong Zhang   }
5083f457be1SHong Zhang 
5094b9ad928SBarry Smith   if (red->pmats) {
5104b9ad928SBarry Smith     ierr = MatDestroyMatrices(1,&red->pmats);CHKERRQ(ierr);
5114b9ad928SBarry Smith   }
5124b9ad928SBarry Smith   ierr = PCDestroy(red->pc);CHKERRQ(ierr);
5134b9ad928SBarry Smith   ierr = PetscFree(red);CHKERRQ(ierr);
5144b9ad928SBarry Smith   PetscFunctionReturn(0);
5154b9ad928SBarry Smith }
5164b9ad928SBarry Smith 
5174b9ad928SBarry Smith #undef __FUNCT__
5184b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant"
5196849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
5204b9ad928SBarry Smith {
5214b9ad928SBarry Smith   PetscFunctionBegin;
5224b9ad928SBarry Smith   PetscFunctionReturn(0);
5234b9ad928SBarry Smith }
5244b9ad928SBarry Smith 
5254b9ad928SBarry Smith EXTERN_C_BEGIN
5264b9ad928SBarry Smith #undef __FUNCT__
5274b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant"
528dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
5294b9ad928SBarry Smith {
5304b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
531dfbe8321SBarry Smith   PetscErrorCode ierr;
5324b9ad928SBarry Smith 
5334b9ad928SBarry Smith   PetscFunctionBegin;
5344b9ad928SBarry Smith   red->scatterin  = in;
5354b9ad928SBarry Smith   red->scatterout = out;
5364b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr);
5374b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr);
5384b9ad928SBarry Smith   PetscFunctionReturn(0);
5394b9ad928SBarry Smith }
5404b9ad928SBarry Smith EXTERN_C_END
5414b9ad928SBarry Smith 
5424b9ad928SBarry Smith #undef __FUNCT__
5434b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter"
5444b9ad928SBarry Smith /*@
5454b9ad928SBarry Smith    PCRedundantSetScatter - Sets the scatter used to copy values into the
5464b9ad928SBarry Smith      redundant local solve and the scatter to move them back into the global
5474b9ad928SBarry Smith      vector.
5484b9ad928SBarry Smith 
5494b9ad928SBarry Smith    Collective on PC
5504b9ad928SBarry Smith 
5514b9ad928SBarry Smith    Input Parameters:
5524b9ad928SBarry Smith +  pc - the preconditioner context
5534b9ad928SBarry Smith .  in - the scatter to move the values in
5544b9ad928SBarry Smith -  out - the scatter to move them out
5554b9ad928SBarry Smith 
5564b9ad928SBarry Smith    Level: advanced
5574b9ad928SBarry Smith 
5584b9ad928SBarry Smith .keywords: PC, redundant solve
5594b9ad928SBarry Smith @*/
560dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
5614b9ad928SBarry Smith {
562dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter);
5634b9ad928SBarry Smith 
5644b9ad928SBarry Smith   PetscFunctionBegin;
5654482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5664482741eSBarry Smith   PetscValidHeaderSpecific(in,VEC_SCATTER_COOKIE,2);
5674482741eSBarry Smith   PetscValidHeaderSpecific(out,VEC_SCATTER_COOKIE,3);
5684b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);CHKERRQ(ierr);
5694b9ad928SBarry Smith   if (f) {
5704b9ad928SBarry Smith     ierr = (*f)(pc,in,out);CHKERRQ(ierr);
5714b9ad928SBarry Smith   }
5724b9ad928SBarry Smith   PetscFunctionReturn(0);
5734b9ad928SBarry Smith }
5744b9ad928SBarry Smith 
5754b9ad928SBarry Smith EXTERN_C_BEGIN
5764b9ad928SBarry Smith #undef __FUNCT__
5774b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC_Redundant"
578dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC_Redundant(PC pc,PC *innerpc)
5794b9ad928SBarry Smith {
5804b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
5814b9ad928SBarry Smith 
5824b9ad928SBarry Smith   PetscFunctionBegin;
5834b9ad928SBarry Smith   *innerpc = red->pc;
5844b9ad928SBarry Smith   PetscFunctionReturn(0);
5854b9ad928SBarry Smith }
5864b9ad928SBarry Smith EXTERN_C_END
5874b9ad928SBarry Smith 
5884b9ad928SBarry Smith #undef __FUNCT__
5894b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC"
5904b9ad928SBarry Smith /*@
5914b9ad928SBarry Smith    PCRedundantGetPC - Gets the sequential PC created by the redundant PC.
5924b9ad928SBarry Smith 
5934b9ad928SBarry Smith    Not Collective
5944b9ad928SBarry Smith 
5954b9ad928SBarry Smith    Input Parameter:
5964b9ad928SBarry Smith .  pc - the preconditioner context
5974b9ad928SBarry Smith 
5984b9ad928SBarry Smith    Output Parameter:
5994b9ad928SBarry Smith .  innerpc - the sequential PC
6004b9ad928SBarry Smith 
6014b9ad928SBarry Smith    Level: advanced
6024b9ad928SBarry Smith 
6034b9ad928SBarry Smith .keywords: PC, redundant solve
6044b9ad928SBarry Smith @*/
605dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC pc,PC *innerpc)
6064b9ad928SBarry Smith {
607dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PC*);
6084b9ad928SBarry Smith 
6094b9ad928SBarry Smith   PetscFunctionBegin;
6104482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6114482741eSBarry Smith   PetscValidPointer(innerpc,2);
6124b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
6134b9ad928SBarry Smith   if (f) {
6144b9ad928SBarry Smith     ierr = (*f)(pc,innerpc);CHKERRQ(ierr);
6154b9ad928SBarry Smith   }
6164b9ad928SBarry Smith   PetscFunctionReturn(0);
6174b9ad928SBarry Smith }
6184b9ad928SBarry Smith 
6194b9ad928SBarry Smith EXTERN_C_BEGIN
6204b9ad928SBarry Smith #undef __FUNCT__
6214b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant"
622dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
6234b9ad928SBarry Smith {
6244b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
6254b9ad928SBarry Smith 
6264b9ad928SBarry Smith   PetscFunctionBegin;
6274b9ad928SBarry Smith   if (mat)  *mat  = red->pmats[0];
6284b9ad928SBarry Smith   if (pmat) *pmat = red->pmats[0];
6294b9ad928SBarry Smith   PetscFunctionReturn(0);
6304b9ad928SBarry Smith }
6314b9ad928SBarry Smith EXTERN_C_END
6324b9ad928SBarry Smith 
6334b9ad928SBarry Smith #undef __FUNCT__
6344b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators"
6354b9ad928SBarry Smith /*@
6364b9ad928SBarry Smith    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
6374b9ad928SBarry Smith 
6384b9ad928SBarry Smith    Not Collective
6394b9ad928SBarry Smith 
6404b9ad928SBarry Smith    Input Parameter:
6414b9ad928SBarry Smith .  pc - the preconditioner context
6424b9ad928SBarry Smith 
6434b9ad928SBarry Smith    Output Parameters:
6444b9ad928SBarry Smith +  mat - the matrix
6454b9ad928SBarry Smith -  pmat - the (possibly different) preconditioner matrix
6464b9ad928SBarry Smith 
6474b9ad928SBarry Smith    Level: advanced
6484b9ad928SBarry Smith 
6494b9ad928SBarry Smith .keywords: PC, redundant solve
6504b9ad928SBarry Smith @*/
651dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
6524b9ad928SBarry Smith {
653dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,Mat*,Mat*);
6544b9ad928SBarry Smith 
6554b9ad928SBarry Smith   PetscFunctionBegin;
6564482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6574482741eSBarry Smith   if (mat)  PetscValidPointer(mat,2);
6584482741eSBarry Smith   if (pmat) PetscValidPointer(pmat,3);
6594b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);CHKERRQ(ierr);
6604b9ad928SBarry Smith   if (f) {
6614b9ad928SBarry Smith     ierr = (*f)(pc,mat,pmat);CHKERRQ(ierr);
6624b9ad928SBarry Smith   }
6634b9ad928SBarry Smith   PetscFunctionReturn(0);
6644b9ad928SBarry Smith }
6654b9ad928SBarry Smith 
6664b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
66737a17b4dSBarry Smith /*MC
66837a17b4dSBarry Smith      PCREDUNDANT - Runs a preconditioner for the entire problem on each processor
66937a17b4dSBarry Smith 
67037a17b4dSBarry Smith 
67137a17b4dSBarry Smith      Options for the redundant preconditioners can be set with -redundant_pc_xxx
67237a17b4dSBarry Smith 
67337a17b4dSBarry Smith    Level: intermediate
67437a17b4dSBarry Smith 
67537a17b4dSBarry Smith 
67637a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
67737a17b4dSBarry Smith            PCRedundantGetPC(), PCRedundantGetOperators()
67837a17b4dSBarry Smith 
67937a17b4dSBarry Smith M*/
68037a17b4dSBarry Smith 
6814b9ad928SBarry Smith EXTERN_C_BEGIN
6824b9ad928SBarry Smith #undef __FUNCT__
6834b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant"
684dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redundant(PC pc)
6854b9ad928SBarry Smith {
686dfbe8321SBarry Smith   PetscErrorCode ierr;
6874b9ad928SBarry Smith   PC_Redundant   *red;
688e060cb09SBarry Smith   const char     *prefix;
6894b9ad928SBarry Smith 
6903f457be1SHong Zhang   PetscInt       nsubcomm,np_subcomm,nleftover,i,j,color;
6913f457be1SHong Zhang   PetscMPIInt    rank,size,subrank,*subsize;
6923f457be1SHong Zhang   MPI_Comm       subcomm;
6933f457be1SHong Zhang   PetscMPIInt    duprank;
6943f457be1SHong Zhang   PetscMPIInt    rank_dup,size_dup;
6953f457be1SHong Zhang   MPI_Comm       dupcomm;
6963f457be1SHong Zhang 
6974b9ad928SBarry Smith   PetscFunctionBegin;
6984b9ad928SBarry Smith   ierr = PetscNew(PC_Redundant,&red);CHKERRQ(ierr);
69952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_Redundant));CHKERRQ(ierr);
7004b9ad928SBarry Smith   red->useparallelmat   = PETSC_TRUE;
7014b9ad928SBarry Smith 
7023f457be1SHong Zhang   ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr);
7033f457be1SHong Zhang   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
7043f457be1SHong Zhang   nsubcomm = size;
7053f457be1SHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-nsubcomm",&nsubcomm,PETSC_NULL);CHKERRQ(ierr);
7063f457be1SHong Zhang   if (nsubcomm > size) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %D cannot be larger than MPI_Comm_size %D",nsubcomm,size);
7073f457be1SHong Zhang 
7083f457be1SHong Zhang   /* get size of each subcommunicators */
7093f457be1SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
7103f457be1SHong Zhang   np_subcomm = size/nsubcomm;
7113f457be1SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
7123f457be1SHong Zhang   for (i=0; i<nsubcomm; i++){
7133f457be1SHong Zhang     subsize[i] = np_subcomm;
7143f457be1SHong Zhang     if (i<nleftover) subsize[i]++;
7153f457be1SHong Zhang   }
7163f457be1SHong Zhang 
7173f457be1SHong Zhang   /* find color for this proc */
7183f457be1SHong Zhang   color   = rank%nsubcomm;
7193f457be1SHong Zhang   subrank = rank/nsubcomm;
7203f457be1SHong Zhang 
7213f457be1SHong Zhang   ierr = MPI_Comm_split(pc->comm,color,subrank,&subcomm);CHKERRQ(ierr);
7223f457be1SHong Zhang   red->subcomm  = subcomm;
7233f457be1SHong Zhang   red->color    = color;
7243f457be1SHong Zhang   red->nsubcomm = nsubcomm;
7253f457be1SHong Zhang 
7263f457be1SHong Zhang   j = 0; duprank = 0;
7273f457be1SHong Zhang   for (i=0; i<nsubcomm; i++){
7283f457be1SHong Zhang     if (j == color){
7293f457be1SHong Zhang       duprank += subrank;
7303f457be1SHong Zhang       break;
7313f457be1SHong Zhang     }
7323f457be1SHong Zhang     duprank += subsize[i]; j++;
7333f457be1SHong Zhang   }
7343f457be1SHong Zhang   /*
7353f457be1SHong Zhang   ierr = PetscSynchronizedPrintf(pc->comm, "[%d] color %d, subrank %d, duprank %d\n",rank,color,subrank,duprank);
7363f457be1SHong Zhang   ierr = PetscSynchronizedFlush(pc->comm);CHKERRQ(ierr);
7373f457be1SHong Zhang   */
7383f457be1SHong Zhang 
7393f457be1SHong Zhang   ierr = MPI_Comm_split(pc->comm,0,duprank,&dupcomm);CHKERRQ(ierr);
7403f457be1SHong Zhang   ierr = MPI_Comm_rank(dupcomm,&rank_dup);CHKERRQ(ierr);
7413f457be1SHong Zhang   ierr = MPI_Comm_size(dupcomm,&size_dup);CHKERRQ(ierr);
7423f457be1SHong Zhang   /*
7433f457be1SHong Zhang   ierr = PetscSynchronizedPrintf(pc->comm, "[%d] duprank %d\n",rank,duprank);
7443f457be1SHong Zhang   ierr = PetscSynchronizedFlush(pc->comm);CHKERRQ(ierr);
7453f457be1SHong Zhang   */
7463f457be1SHong Zhang   red->dupcomm = dupcomm;
7473f457be1SHong Zhang   ierr = PetscFree(subsize);CHKERRQ(ierr);
7483f457be1SHong Zhang 
7494b9ad928SBarry Smith   /* create the sequential PC that each processor has copy of */
7503f457be1SHong Zhang   ierr = PCCreate(subcomm,&red->pc);CHKERRQ(ierr);
7514b9ad928SBarry Smith   ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr);
7524b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
7534b9ad928SBarry Smith   ierr = PCSetOptionsPrefix(red->pc,prefix);CHKERRQ(ierr);
7544b9ad928SBarry Smith   ierr = PCAppendOptionsPrefix(red->pc,"redundant_");CHKERRQ(ierr);
7554b9ad928SBarry Smith 
7564b9ad928SBarry Smith   pc->ops->apply             = PCApply_Redundant;
7574b9ad928SBarry Smith   pc->ops->applytranspose    = 0;
7584b9ad928SBarry Smith   pc->ops->setup             = PCSetUp_Redundant;
7594b9ad928SBarry Smith   pc->ops->destroy           = PCDestroy_Redundant;
7604b9ad928SBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_Redundant;
7614b9ad928SBarry Smith   pc->ops->setuponblocks     = 0;
7624b9ad928SBarry Smith   pc->ops->view              = PCView_Redundant;
7634b9ad928SBarry Smith   pc->ops->applyrichardson   = 0;
7644b9ad928SBarry Smith 
7654b9ad928SBarry Smith   pc->data     = (void*)red;
7664b9ad928SBarry Smith 
7674b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
7684b9ad928SBarry Smith                     PCRedundantSetScatter_Redundant);CHKERRQ(ierr);
7694b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant",
7704b9ad928SBarry Smith                     PCRedundantGetPC_Redundant);CHKERRQ(ierr);
7714b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
7724b9ad928SBarry Smith                     PCRedundantGetOperators_Redundant);CHKERRQ(ierr);
7734b9ad928SBarry Smith 
7744b9ad928SBarry Smith   PetscFunctionReturn(0);
7754b9ad928SBarry Smith }
7764b9ad928SBarry Smith EXTERN_C_END
777