xref: /petsc/src/ksp/pc/impls/redistribute/redistribute.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
1df826632SBarry Smith 
2df826632SBarry Smith /*
3df826632SBarry Smith   This file defines a "solve the problem redistributely on each subgroup of processor" preconditioner.
4df826632SBarry Smith */
5b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscksp.h" I*/
6c6db04a5SJed Brown #include <petscksp.h>
7df826632SBarry Smith 
8df826632SBarry Smith typedef struct {
9df826632SBarry Smith   KSP         ksp;
10df826632SBarry Smith   Vec         x,b;
11df826632SBarry Smith   VecScatter  scatter;
12911f9fe8SBarry Smith   IS          is;
13181dd334SBarry Smith   PetscInt    dcnt,*drows;    /* these are the local rows that have only diagonal entry */
14181dd334SBarry Smith   PetscScalar *diag;
15181dd334SBarry Smith   Vec         work;
16df826632SBarry Smith } PC_Redistribute;
17df826632SBarry Smith 
18df826632SBarry Smith #undef __FUNCT__
19df826632SBarry Smith #define __FUNCT__ "PCView_Redistribute"
20df826632SBarry Smith static PetscErrorCode PCView_Redistribute(PC pc,PetscViewer viewer)
21df826632SBarry Smith {
22df826632SBarry Smith   PC_Redistribute *red = (PC_Redistribute*)pc->data;
23df826632SBarry Smith   PetscErrorCode  ierr;
24ace3abfcSBarry Smith   PetscBool       iascii,isstring;
25181dd334SBarry Smith   PetscInt        ncnt,N;
26df826632SBarry Smith 
27df826632SBarry Smith   PetscFunctionBegin;
28251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
30df826632SBarry Smith   if (iascii) {
31181dd334SBarry Smith     ierr = MPI_Allreduce(&red->dcnt,&ncnt,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
32181dd334SBarry Smith     ierr = MatGetSize(pc->pmat,&N,PETSC_NULL);CHKERRQ(ierr);
33572f72c7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"    Number rows eliminated %D Percentage rows eliminated %g\n",ncnt,100.0*((PetscReal)ncnt)/((PetscReal)N));CHKERRQ(ierr);
34df826632SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Redistribute preconditioner: \n");CHKERRQ(ierr);
35df826632SBarry Smith     ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr);
36df826632SBarry Smith   } else if (isstring) {
37df826632SBarry Smith     ierr = PetscViewerStringSPrintf(viewer," Redistribute preconditioner");CHKERRQ(ierr);
38df826632SBarry Smith     ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr);
3911aeaf0aSBarry Smith   }
40df826632SBarry Smith   PetscFunctionReturn(0);
41df826632SBarry Smith }
42df826632SBarry Smith 
43df826632SBarry Smith #undef __FUNCT__
44df826632SBarry Smith #define __FUNCT__ "PCSetUp_Redistribute"
45df826632SBarry Smith static PetscErrorCode PCSetUp_Redistribute(PC pc)
46df826632SBarry Smith {
47df826632SBarry Smith   PC_Redistribute   *red = (PC_Redistribute*)pc->data;
48df826632SBarry Smith   PetscErrorCode    ierr;
49df826632SBarry Smith   MPI_Comm          comm;
50181dd334SBarry Smith   PetscInt          rstart,rend,i,nz,cnt,*rows,ncnt,dcnt,*drows;
5126283091SBarry Smith   PetscLayout       map,nmap;
525a586d82SBarry Smith   PetscMPIInt       size,imdex,tag,n;
53e8dd6687SHong Zhang   PetscInt          *source = PETSC_NULL;
54df826632SBarry Smith   PetscMPIInt       *nprocs = PETSC_NULL,nrecvs;
55df826632SBarry Smith   PetscInt          j,nsends;
56df826632SBarry Smith   PetscInt          *owner = PETSC_NULL,*starts = PETSC_NULL,count,slen;
57e8dd6687SHong Zhang   PetscInt          *rvalues,*svalues,recvtotal;
58df826632SBarry Smith   PetscMPIInt       *onodes1,*olengths1;
59df826632SBarry Smith   MPI_Request       *send_waits = PETSC_NULL,*recv_waits = PETSC_NULL;
60df826632SBarry Smith   MPI_Status        recv_status,*send_status;
61181dd334SBarry Smith   Vec               tvec,diag;
62ca320bd4SBarry Smith   Mat               tmat;
63181dd334SBarry Smith   const PetscScalar *d;
64df826632SBarry Smith 
65df826632SBarry Smith   PetscFunctionBegin;
66dc9360f3SBarry Smith   if (pc->setupcalled) {
67dc9360f3SBarry Smith     ierr = KSPGetOperators(red->ksp,PETSC_NULL,&tmat,PETSC_NULL);CHKERRQ(ierr);
68dc9360f3SBarry Smith     ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_REUSE_MATRIX,&tmat);CHKERRQ(ierr);
691d805cfdSBarry Smith     ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
70dc9360f3SBarry Smith   } else {
71b862ddfaSBarry Smith     PetscInt NN;
72b862ddfaSBarry Smith 
73df826632SBarry Smith     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
74df826632SBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
75361c1e09SMatthew Knepley     ierr = PetscObjectGetNewTag((PetscObject)pc,&tag);CHKERRQ(ierr);
76df826632SBarry Smith 
77ca320bd4SBarry Smith     /* count non-diagonal rows on process */
78df826632SBarry Smith     ierr = MatGetOwnershipRange(pc->mat,&rstart,&rend);CHKERRQ(ierr);
79ca320bd4SBarry Smith     cnt  = 0;
80df826632SBarry Smith     for (i=rstart; i<rend; i++) {
81df826632SBarry Smith       ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
82df826632SBarry Smith       if (nz > 1) cnt++;
83911f9fe8SBarry Smith       ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
84df826632SBarry Smith     }
85df826632SBarry Smith     ierr = PetscMalloc(cnt*sizeof(PetscInt),&rows);CHKERRQ(ierr);
86181dd334SBarry Smith     ierr = PetscMalloc((rend - rstart - cnt)*sizeof(PetscInt),&drows);CHKERRQ(ierr);
87ca320bd4SBarry Smith 
88ca320bd4SBarry Smith     /* list non-diagonal rows on process */
89181dd334SBarry Smith     cnt = 0; dcnt = 0;
90df826632SBarry Smith     for (i=rstart; i<rend; i++) {
91df826632SBarry Smith       ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
92df826632SBarry Smith       if (nz > 1) rows[cnt++] = i;
93181dd334SBarry Smith       else drows[dcnt++] = i - rstart;
94911f9fe8SBarry Smith       ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
95df826632SBarry Smith     }
96ca320bd4SBarry Smith 
9726283091SBarry Smith     /* create PetscLayout for non-diagonal rows on each process */
9826283091SBarry Smith     ierr   = PetscLayoutCreate(comm,&map);CHKERRQ(ierr);
9926283091SBarry Smith     ierr   = PetscLayoutSetLocalSize(map,cnt);CHKERRQ(ierr);
10026283091SBarry Smith     ierr   = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
10126283091SBarry Smith     ierr   = PetscLayoutSetUp(map);CHKERRQ(ierr);
102df826632SBarry Smith     rstart = map->rstart;
103df826632SBarry Smith     rend   = map->rend;
104df826632SBarry Smith 
10526283091SBarry Smith     /* create PetscLayout for load-balanced non-diagonal rows on each process */
10626283091SBarry Smith     ierr = PetscLayoutCreate(comm,&nmap);CHKERRQ(ierr);
107df826632SBarry Smith     ierr = MPI_Allreduce(&cnt,&ncnt,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
10826283091SBarry Smith     ierr = PetscLayoutSetSize(nmap,ncnt);CHKERRQ(ierr);
10926283091SBarry Smith     ierr = PetscLayoutSetBlockSize(nmap,1);CHKERRQ(ierr);
11026283091SBarry Smith     ierr = PetscLayoutSetUp(nmap);CHKERRQ(ierr);
111df826632SBarry Smith 
112b862ddfaSBarry Smith     ierr = MatGetSize(pc->pmat,&NN,PETSC_NULL);CHKERRQ(ierr);
1133bf036e2SBarry Smith     ierr = PetscInfo2(pc,"Number of diagonal rows eliminated %d, percentage eliminated %g\n",NN-ncnt,((PetscReal)(NN-ncnt))/((PetscReal)(NN)));CHKERRQ(ierr);
114ca320bd4SBarry Smith     /*
115ca320bd4SBarry Smith         this code is taken from VecScatterCreate_PtoS()
116ca320bd4SBarry Smith         Determines what rows need to be moved where to
117ca320bd4SBarry Smith         load balance the non-diagonal rows
118ca320bd4SBarry Smith     */
119df826632SBarry Smith     /*  count number of contributors to each processor */
120df826632SBarry Smith     ierr   = PetscMalloc2(size,PetscMPIInt,&nprocs,cnt,PetscInt,&owner);CHKERRQ(ierr);
121df826632SBarry Smith     ierr   = PetscMemzero(nprocs,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
122df826632SBarry Smith     j      = 0;
123df826632SBarry Smith     nsends = 0;
124df826632SBarry Smith     for (i=rstart; i<rend; i++) {
125df826632SBarry Smith       if (i < nmap->range[j]) j = 0;
126df826632SBarry Smith       for (; j<size; j++) {
127df826632SBarry Smith         if (i < nmap->range[j+1]) {
128df826632SBarry Smith           if (!nprocs[j]++) nsends++;
129ca320bd4SBarry Smith           owner[i-rstart] = j;
130df826632SBarry Smith           break;
131df826632SBarry Smith         }
132df826632SBarry Smith       }
133df826632SBarry Smith     }
134df826632SBarry Smith     /* inform other processors of number of messages and max length*/
135df826632SBarry Smith     ierr      = PetscGatherNumberOfMessages(comm,PETSC_NULL,nprocs,&nrecvs);CHKERRQ(ierr);
136df826632SBarry Smith     ierr      = PetscGatherMessageLengths(comm,nsends,nrecvs,nprocs,&onodes1,&olengths1);CHKERRQ(ierr);
137df826632SBarry Smith     ierr      = PetscSortMPIIntWithArray(nrecvs,onodes1,olengths1);CHKERRQ(ierr);
138df826632SBarry Smith     recvtotal = 0; for (i=0; i<nrecvs; i++) recvtotal += olengths1[i];
139df826632SBarry Smith 
140df826632SBarry Smith     /* post receives:  rvalues - rows I will own; count - nu */
141df826632SBarry Smith     ierr  = PetscMalloc3(recvtotal,PetscInt,&rvalues,nrecvs,PetscInt,&source,nrecvs,MPI_Request,&recv_waits);CHKERRQ(ierr);
142df826632SBarry Smith     count = 0;
143df826632SBarry Smith     for (i=0; i<nrecvs; i++) {
144df826632SBarry Smith       ierr   = MPI_Irecv((rvalues+count),olengths1[i],MPIU_INT,onodes1[i],tag,comm,recv_waits+i);CHKERRQ(ierr);
145df826632SBarry Smith       count += olengths1[i];
146df826632SBarry Smith     }
147df826632SBarry Smith 
148df826632SBarry Smith     /* do sends:
149df826632SBarry Smith        1) starts[i] gives the starting index in svalues for stuff going to
150df826632SBarry Smith        the ith processor
151df826632SBarry Smith     */
152911f9fe8SBarry Smith     ierr      = PetscMalloc3(cnt,PetscInt,&svalues,nsends,MPI_Request,&send_waits,size,PetscInt,&starts);CHKERRQ(ierr);
153df826632SBarry Smith     starts[0] = 0;
154*2fa5cd67SKarl Rupp     for (i=1; i<size; i++) starts[i] = starts[i-1] + nprocs[i-1];
155*2fa5cd67SKarl Rupp     for (i=0; i<cnt; i++)  svalues[starts[owner[i]]++] = rows[i];
156181dd334SBarry Smith     for (i=0; i<cnt; i++)  rows[i] = rows[i] - rstart;
157181dd334SBarry Smith     red->drows = drows;
158181dd334SBarry Smith     red->dcnt  = dcnt;
159ca320bd4SBarry Smith     ierr       = PetscFree(rows);CHKERRQ(ierr);
160181dd334SBarry Smith 
161df826632SBarry Smith     starts[0] = 0;
162*2fa5cd67SKarl Rupp     for (i=1; i<size; i++) starts[i] = starts[i-1] + nprocs[i-1];
163df826632SBarry Smith     count = 0;
164df826632SBarry Smith     for (i=0; i<size; i++) {
165df826632SBarry Smith       if (nprocs[i]) {
166df826632SBarry Smith         ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
167df826632SBarry Smith       }
168df826632SBarry Smith     }
169df826632SBarry Smith 
170df826632SBarry Smith     /*  wait on receives */
171df826632SBarry Smith     count = nrecvs;
172df826632SBarry Smith     slen  = 0;
173df826632SBarry Smith     while (count) {
174df826632SBarry Smith       ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
175df826632SBarry Smith       /* unpack receives into our local space */
176df826632SBarry Smith       ierr  = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
177df826632SBarry Smith       slen += n;
178df826632SBarry Smith       count--;
179df826632SBarry Smith     }
180e32f2f54SBarry Smith     if (slen != recvtotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Total message lengths %D not expected %D",slen,recvtotal);
181df826632SBarry Smith 
18270b3c8c7SBarry Smith     ierr = ISCreateGeneral(comm,slen,rvalues,PETSC_COPY_VALUES,&red->is);CHKERRQ(ierr);
183911f9fe8SBarry Smith 
184ca320bd4SBarry Smith     /* free up all work space */
185df826632SBarry Smith     ierr = PetscFree(olengths1);CHKERRQ(ierr);
186df826632SBarry Smith     ierr = PetscFree(onodes1);CHKERRQ(ierr);
187df826632SBarry Smith     ierr = PetscFree3(rvalues,source,recv_waits);CHKERRQ(ierr);
188ca320bd4SBarry Smith     ierr = PetscFree2(nprocs,owner);CHKERRQ(ierr);
189ca320bd4SBarry Smith     if (nsends) {   /* wait on sends */
190df826632SBarry Smith       ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
191df826632SBarry Smith       ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
192df826632SBarry Smith       ierr = PetscFree(send_status);CHKERRQ(ierr);
193df826632SBarry Smith     }
194df826632SBarry Smith     ierr = PetscFree3(svalues,send_waits,starts);CHKERRQ(ierr);
1956bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1966bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&nmap);CHKERRQ(ierr);
197df826632SBarry Smith 
198ca320bd4SBarry Smith     ierr = VecCreateMPI(comm,slen,PETSC_DETERMINE,&red->b);CHKERRQ(ierr);
199ca320bd4SBarry Smith     ierr = VecDuplicate(red->b,&red->x);CHKERRQ(ierr);
200ca320bd4SBarry Smith     ierr = MatGetVecs(pc->pmat,&tvec,PETSC_NULL);CHKERRQ(ierr);
201ca320bd4SBarry Smith     ierr = VecScatterCreate(tvec,red->is,red->b,PETSC_NULL,&red->scatter);CHKERRQ(ierr);
2026bf464f9SBarry Smith     ierr = VecDestroy(&tvec);CHKERRQ(ierr);
203ca320bd4SBarry Smith     ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_INITIAL_MATRIX,&tmat);CHKERRQ(ierr);
204ca320bd4SBarry Smith     ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2056bf464f9SBarry Smith     ierr = MatDestroy(&tmat);CHKERRQ(ierr);
2061d805cfdSBarry Smith   }
207181dd334SBarry Smith 
208181dd334SBarry Smith   /* get diagonal portion of matrix */
209181dd334SBarry Smith   ierr = PetscMalloc(red->dcnt*sizeof(PetscScalar),&red->diag);CHKERRQ(ierr);
210181dd334SBarry Smith   ierr = MatGetVecs(pc->pmat,&diag,PETSC_NULL);CHKERRQ(ierr);
211181dd334SBarry Smith   ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr);
2123649974fSBarry Smith   ierr = VecGetArrayRead(diag,&d);CHKERRQ(ierr);
213*2fa5cd67SKarl Rupp   for (i=0; i<red->dcnt; i++) red->diag[i] = 1.0/d[red->drows[i]];
2143649974fSBarry Smith   ierr = VecRestoreArrayRead(diag,&d);CHKERRQ(ierr);
2156bf464f9SBarry Smith   ierr = VecDestroy(&diag);CHKERRQ(ierr);
216df826632SBarry Smith   ierr = KSPSetUp(red->ksp);CHKERRQ(ierr);
217df826632SBarry Smith   PetscFunctionReturn(0);
218df826632SBarry Smith }
219df826632SBarry Smith 
220df826632SBarry Smith #undef __FUNCT__
221df826632SBarry Smith #define __FUNCT__ "PCApply_Redistribute"
222df826632SBarry Smith static PetscErrorCode PCApply_Redistribute(PC pc,Vec b,Vec x)
223df826632SBarry Smith {
224df826632SBarry Smith   PC_Redistribute   *red = (PC_Redistribute*)pc->data;
225df826632SBarry Smith   PetscErrorCode    ierr;
226181dd334SBarry Smith   PetscInt          dcnt   = red->dcnt,i;
227181dd334SBarry Smith   const PetscInt    *drows = red->drows;
228181dd334SBarry Smith   PetscScalar       *xwork;
229181dd334SBarry Smith   const PetscScalar *bwork,*diag = red->diag;
230df826632SBarry Smith 
231df826632SBarry Smith   PetscFunctionBegin;
232181dd334SBarry Smith   if (!red->work) {
233181dd334SBarry Smith     ierr = VecDuplicate(b,&red->work);CHKERRQ(ierr);
234181dd334SBarry Smith   }
235181dd334SBarry Smith   /* compute the rows of solution that have diagonal entries only */
236181dd334SBarry Smith   ierr = VecSet(x,0.0);CHKERRQ(ierr);         /* x = diag(A)^{-1} b */
237181dd334SBarry Smith   ierr = VecGetArray(x,&xwork);CHKERRQ(ierr);
2383649974fSBarry Smith   ierr = VecGetArrayRead(b,&bwork);CHKERRQ(ierr);
239*2fa5cd67SKarl Rupp   for (i=0; i<dcnt; i++) xwork[drows[i]] = diag[i]*bwork[drows[i]];
240181dd334SBarry Smith   ierr = PetscLogFlops(dcnt);CHKERRQ(ierr);
241181dd334SBarry Smith   ierr = VecRestoreArray(red->work,&xwork);CHKERRQ(ierr);
2423649974fSBarry Smith   ierr = VecRestoreArrayRead(b,&bwork);CHKERRQ(ierr);
243181dd334SBarry Smith   /* update the right hand side for the reduced system with diagonal rows (and corresponding columns) removed */
244181dd334SBarry Smith   ierr = MatMult(pc->pmat,x,red->work);CHKERRQ(ierr);
245181dd334SBarry Smith   ierr = VecAYPX(red->work,-1.0,b);CHKERRQ(ierr);   /* red->work = b - A x */
246181dd334SBarry Smith 
247181dd334SBarry Smith   ierr = VecScatterBegin(red->scatter,red->work,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
248181dd334SBarry Smith   ierr = VecScatterEnd(red->scatter,red->work,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
249df826632SBarry Smith   ierr = KSPSolve(red->ksp,red->b,red->x);CHKERRQ(ierr);
250df826632SBarry Smith   ierr = VecScatterBegin(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
251df826632SBarry Smith   ierr = VecScatterEnd(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
252df826632SBarry Smith   PetscFunctionReturn(0);
253df826632SBarry Smith }
254df826632SBarry Smith 
255df826632SBarry Smith #undef __FUNCT__
256df826632SBarry Smith #define __FUNCT__ "PCDestroy_Redistribute"
257df826632SBarry Smith static PetscErrorCode PCDestroy_Redistribute(PC pc)
258df826632SBarry Smith {
259df826632SBarry Smith   PC_Redistribute *red = (PC_Redistribute*)pc->data;
260df826632SBarry Smith   PetscErrorCode  ierr;
261df826632SBarry Smith 
262df826632SBarry Smith   PetscFunctionBegin;
263fcfd50ebSBarry Smith   ierr = VecScatterDestroy(&red->scatter);CHKERRQ(ierr);
264fcfd50ebSBarry Smith   ierr = ISDestroy(&red->is);CHKERRQ(ierr);
265fcfd50ebSBarry Smith   ierr = VecDestroy(&red->b);CHKERRQ(ierr);
266fcfd50ebSBarry Smith   ierr = VecDestroy(&red->x);CHKERRQ(ierr);
267fcfd50ebSBarry Smith   ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr);
268fcfd50ebSBarry Smith   ierr = VecDestroy(&red->work);CHKERRQ(ierr);
2693bf036e2SBarry Smith   ierr = PetscFree(red->drows);CHKERRQ(ierr);
2703bf036e2SBarry Smith   ierr = PetscFree(red->diag);CHKERRQ(ierr);
271c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
272df826632SBarry Smith   PetscFunctionReturn(0);
273df826632SBarry Smith }
274df826632SBarry Smith 
275df826632SBarry Smith #undef __FUNCT__
276df826632SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redistribute"
277df826632SBarry Smith static PetscErrorCode PCSetFromOptions_Redistribute(PC pc)
278df826632SBarry Smith {
279df826632SBarry Smith   PetscErrorCode  ierr;
280df826632SBarry Smith   PC_Redistribute *red = (PC_Redistribute*)pc->data;
281df826632SBarry Smith 
282df826632SBarry Smith   PetscFunctionBegin;
283df826632SBarry Smith   ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr);
284df826632SBarry Smith   PetscFunctionReturn(0);
285df826632SBarry Smith }
286df826632SBarry Smith 
2875e7ef714SBarry Smith #undef __FUNCT__
2885e7ef714SBarry Smith #define __FUNCT__ "PCRedistributeGetKSP"
2895e7ef714SBarry Smith /*@
2905e7ef714SBarry Smith    PCRedistributeGetKSP - Gets the KSP created by the PCREDISTRIBUTE
2915e7ef714SBarry Smith 
2925e7ef714SBarry Smith    Not Collective
2935e7ef714SBarry Smith 
2945e7ef714SBarry Smith    Input Parameter:
2955e7ef714SBarry Smith .  pc - the preconditioner context
2965e7ef714SBarry Smith 
2975e7ef714SBarry Smith    Output Parameter:
2985e7ef714SBarry Smith .  innerksp - the inner KSP
2995e7ef714SBarry Smith 
3005e7ef714SBarry Smith    Level: advanced
3015e7ef714SBarry Smith 
3025e7ef714SBarry Smith .keywords: PC, redistribute solve
3035e7ef714SBarry Smith @*/
3045e7ef714SBarry Smith PetscErrorCode  PCRedistributeGetKSP(PC pc,KSP *innerksp)
3055e7ef714SBarry Smith {
3065e7ef714SBarry Smith   PC_Redistribute *red = (PC_Redistribute*)pc->data;
3075e7ef714SBarry Smith 
3085e7ef714SBarry Smith   PetscFunctionBegin;
3095e7ef714SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3105e7ef714SBarry Smith   PetscValidPointer(innerksp,2);
3115e7ef714SBarry Smith   *innerksp = red->ksp;
3125e7ef714SBarry Smith   PetscFunctionReturn(0);
3135e7ef714SBarry Smith }
3145e7ef714SBarry Smith 
315df826632SBarry Smith /* -------------------------------------------------------------------------------------*/
316df826632SBarry Smith /*MC
317181dd334SBarry Smith      PCREDISTRIBUTE - Redistributes a matrix for load balancing, removing the rows that only have a diagonal entry and then applys a KSP to that new matrix
318df826632SBarry Smith 
319df826632SBarry Smith      Options for the redistribute preconditioners can be set with -redistribute_ksp_xxx <values> and -redistribute_pc_xxx <values>
320df826632SBarry Smith 
321181dd334SBarry Smith      Notes:  Usually run this with -ksp_type preonly
322181dd334SBarry Smith 
323181dd334SBarry Smith      If you have used MatZeroRows() to eliminate (for example, Dirichlet) boundary conditions for a symmetric problem then you can use, for example, -ksp_type preonly
324181dd334SBarry Smith      -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc to take advantage of the symmetry.
325181dd334SBarry Smith 
3262dfef595SBarry Smith      This does NOT call a partitioner to reorder rows to lower communication; the ordering of the rows in the original matrix and redistributed matrix is the same.
3272dfef595SBarry Smith 
3282dfef595SBarry Smith      Developer Notes: Should add an option to this preconditioner to use a partitioner to redistribute the rows to lower communication.
3292dfef595SBarry Smith 
330df826632SBarry Smith    Level: intermediate
331df826632SBarry Smith 
3325e7ef714SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedistributeGetKSP()
333df826632SBarry Smith M*/
334df826632SBarry Smith 
335df826632SBarry Smith EXTERN_C_BEGIN
336df826632SBarry Smith #undef __FUNCT__
337df826632SBarry Smith #define __FUNCT__ "PCCreate_Redistribute"
3387087cfbeSBarry Smith PetscErrorCode  PCCreate_Redistribute(PC pc)
339df826632SBarry Smith {
340df826632SBarry Smith   PetscErrorCode  ierr;
341df826632SBarry Smith   PC_Redistribute *red;
342911f9fe8SBarry Smith   const char      *prefix;
343df826632SBarry Smith 
344df826632SBarry Smith   PetscFunctionBegin;
345df826632SBarry Smith   ierr     = PetscNewLog(pc,PC_Redistribute,&red);CHKERRQ(ierr);
346df826632SBarry Smith   pc->data = (void*)red;
347df826632SBarry Smith 
348df826632SBarry Smith   pc->ops->apply          = PCApply_Redistribute;
349df826632SBarry Smith   pc->ops->applytranspose = 0;
350df826632SBarry Smith   pc->ops->setup          = PCSetUp_Redistribute;
351df826632SBarry Smith   pc->ops->destroy        = PCDestroy_Redistribute;
352df826632SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Redistribute;
353df826632SBarry Smith   pc->ops->view           = PCView_Redistribute;
354911f9fe8SBarry Smith 
355911f9fe8SBarry Smith   ierr = KSPCreate(((PetscObject)pc)->comm,&red->ksp);CHKERRQ(ierr);
356911f9fe8SBarry Smith   ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
357911f9fe8SBarry Smith   ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr);
358911f9fe8SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
359911f9fe8SBarry Smith   ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
360911f9fe8SBarry Smith   ierr = KSPAppendOptionsPrefix(red->ksp,"redistribute_");CHKERRQ(ierr);
361df826632SBarry Smith   PetscFunctionReturn(0);
362df826632SBarry Smith }
363df826632SBarry Smith EXTERN_C_END
364