xref: /petsc/src/ksp/pc/impls/redistribute/redistribute.c (revision f89e2e87fbc2bf363b9a1b6d11cc8db8663cdbc0)
1 #define PETSCKSP_DLL
2 
3 /*
4   This file defines a "solve the problem redistributely on each subgroup of processor" preconditioner.
5 */
6 #include "private/pcimpl.h"     /*I "petscpc.h" I*/
7 #include "petscksp.h"
8 
9 typedef struct {
10   KSP          ksp;
11   Vec          x,b;
12   VecScatter   scatter;
13   IS           is;
14   PetscInt     dcnt,*drows;   /* these are the local rows that have only diagonal entry */
15   PetscScalar  *diag;
16   Vec          work;
17 } PC_Redistribute;
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "PCView_Redistribute"
21 static PetscErrorCode PCView_Redistribute(PC pc,PetscViewer viewer)
22 {
23   PC_Redistribute *red = (PC_Redistribute*)pc->data;
24   PetscErrorCode  ierr;
25   PetscTruth      iascii,isstring;
26   PetscInt        ncnt,N;
27 
28   PetscFunctionBegin;
29   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
30   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
31   if (iascii) {
32     ierr = MPI_Allreduce(&red->dcnt,&ncnt,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
33     ierr = MatGetSize(pc->pmat,&N,PETSC_NULL);CHKERRQ(ierr);
34     ierr = PetscViewerASCIIPrintf(viewer,"    Number rows eliminated %D Percentage rows eliminated %g\n",ncnt,100.0*((PetscReal)ncnt)/((PetscReal)N));CHKERRQ(ierr);
35     ierr = PetscViewerASCIIPrintf(viewer,"  Redistribute preconditioner: \n");CHKERRQ(ierr);
36     ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr);
37   } else if (isstring) {
38     ierr = PetscViewerStringSPrintf(viewer," Redistribute preconditioner");CHKERRQ(ierr);
39     ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr);
40   } else {
41     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PC redistribute",((PetscObject)viewer)->type_name);
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "PCSetUp_Redistribute"
48 static PetscErrorCode PCSetUp_Redistribute(PC pc)
49 {
50   PC_Redistribute   *red = (PC_Redistribute*)pc->data;
51   PetscErrorCode    ierr;
52   MPI_Comm          comm;
53   PetscInt          rstart,rend,i,nz,cnt,*rows,ncnt,dcnt,*drows;
54   PetscLayout       map,nmap;
55   PetscMPIInt       size,rank,imdex,tag,n;
56   PetscInt          *source = PETSC_NULL;
57   PetscMPIInt       *nprocs = PETSC_NULL,nrecvs;
58   PetscInt          j,nsends;
59   PetscInt          *owner = PETSC_NULL,*starts = PETSC_NULL,count,slen;
60   PetscInt          *rvalues,*svalues,recvtotal;
61   PetscMPIInt       *onodes1,*olengths1;
62   MPI_Request       *send_waits = PETSC_NULL,*recv_waits = PETSC_NULL;
63   MPI_Status        recv_status,*send_status;
64   Vec               tvec,diag;
65   Mat               tmat;
66   const PetscScalar *d;
67 
68   PetscFunctionBegin;
69   if (pc->setupcalled) {
70     ierr = KSPGetOperators(red->ksp,PETSC_NULL,&tmat,PETSC_NULL);CHKERRQ(ierr);
71     ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_REUSE_MATRIX,&tmat);CHKERRQ(ierr);
72     ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
73   } else {
74     PetscInt NN;
75 
76     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
77     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
78     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
79     ierr = PetscObjectGetNewTag((PetscObject)pc,&tag);CHKERRQ(ierr);
80 
81     /* count non-diagonal rows on process */
82     ierr = MatGetOwnershipRange(pc->mat,&rstart,&rend);CHKERRQ(ierr);
83     cnt  = 0;
84     for (i=rstart; i<rend; i++) {
85       ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
86       if (nz > 1) cnt++;
87       ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
88     }
89     ierr = PetscMalloc(cnt*sizeof(PetscInt),&rows);CHKERRQ(ierr);
90     ierr = PetscMalloc((rend - rstart - cnt)*sizeof(PetscInt),&drows);CHKERRQ(ierr);
91 
92     /* list non-diagonal rows on process */
93     cnt  = 0; dcnt = 0;
94     for (i=rstart; i<rend; i++) {
95       ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
96       if (nz > 1) rows[cnt++] = i;
97       else drows[dcnt++] = i - rstart;
98       ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
99     }
100 
101     /* create PetscLayout for non-diagonal rows on each process */
102     ierr = PetscLayoutCreate(comm,&map);CHKERRQ(ierr);
103     ierr = PetscLayoutSetLocalSize(map,cnt);CHKERRQ(ierr);
104     ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
105     ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
106     rstart = map->rstart;
107     rend   = map->rend;
108 
109     /* create PetscLayout for load-balanced non-diagonal rows on each process */
110     ierr = PetscLayoutCreate(comm,&nmap);CHKERRQ(ierr);
111     ierr = MPI_Allreduce(&cnt,&ncnt,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
112     ierr = PetscLayoutSetSize(nmap,ncnt);CHKERRQ(ierr);
113     ierr = PetscLayoutSetBlockSize(nmap,1);CHKERRQ(ierr);
114     ierr = PetscLayoutSetUp(nmap);CHKERRQ(ierr);
115 
116     ierr = MatGetSize(pc->pmat,&NN,PETSC_NULL);CHKERRQ(ierr);
117     ierr = PetscInfo2(pc,"Number of diagonal rows eliminated %d, percentage eliminated %g\n",NN-ncnt,((PetscReal)(NN-ncnt))/((PetscReal)(NN)));
118     /*
119 	this code is taken from VecScatterCreate_PtoS()
120 	Determines what rows need to be moved where to
121 	load balance the non-diagonal rows
122     */
123     /*  count number of contributors to each processor */
124     ierr = PetscMalloc2(size,PetscMPIInt,&nprocs,cnt,PetscInt,&owner);CHKERRQ(ierr);
125     ierr = PetscMemzero(nprocs,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
126     j      = 0;
127     nsends = 0;
128     for (i=rstart; i<rend; i++) {
129       if (i < nmap->range[j]) j = 0;
130       for (; j<size; j++) {
131 	if (i < nmap->range[j+1]) {
132 	  if (!nprocs[j]++) nsends++;
133 	  owner[i-rstart] = j;
134 	  break;
135 	}
136       }
137     }
138     /* inform other processors of number of messages and max length*/
139     ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,nprocs,&nrecvs);CHKERRQ(ierr);
140     ierr = PetscGatherMessageLengths(comm,nsends,nrecvs,nprocs,&onodes1,&olengths1);CHKERRQ(ierr);
141     ierr = PetscSortMPIIntWithArray(nrecvs,onodes1,olengths1);CHKERRQ(ierr);
142     recvtotal = 0; for (i=0; i<nrecvs; i++) recvtotal += olengths1[i];
143 
144     /* post receives:  rvalues - rows I will own; count - nu */
145     ierr = PetscMalloc3(recvtotal,PetscInt,&rvalues,nrecvs,PetscInt,&source,nrecvs,MPI_Request,&recv_waits);CHKERRQ(ierr);
146     count  = 0;
147     for (i=0; i<nrecvs; i++) {
148       ierr  = MPI_Irecv((rvalues+count),olengths1[i],MPIU_INT,onodes1[i],tag,comm,recv_waits+i);CHKERRQ(ierr);
149       count += olengths1[i];
150     }
151 
152     /* do sends:
153        1) starts[i] gives the starting index in svalues for stuff going to
154        the ith processor
155     */
156     ierr = PetscMalloc3(cnt,PetscInt,&svalues,nsends,MPI_Request,&send_waits,size,PetscInt,&starts);CHKERRQ(ierr);
157     starts[0]  = 0;
158     for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
159     for (i=0; i<cnt; i++) {
160       svalues[starts[owner[i]]++] = rows[i];
161     }
162     for (i=0; i<cnt; i++) rows[i] = rows[i] - rstart;
163     red->drows = drows;
164     red->dcnt  = dcnt;
165     ierr = PetscFree(rows);CHKERRQ(ierr);
166 
167     starts[0] = 0;
168     for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
169     count = 0;
170     for (i=0; i<size; i++) {
171       if (nprocs[i]) {
172 	ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
173       }
174     }
175 
176     /*  wait on receives */
177     count  = nrecvs;
178     slen   = 0;
179     while (count) {
180       ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
181       /* unpack receives into our local space */
182       ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
183       slen += n;
184       count--;
185     }
186     if (slen != recvtotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Total message lengths %D not expected %D",slen,recvtotal);
187 
188     ierr = ISCreateGeneral(comm,slen,rvalues,&red->is);CHKERRQ(ierr);
189 
190     /* free up all work space */
191     ierr = PetscFree(olengths1);CHKERRQ(ierr);
192     ierr = PetscFree(onodes1);CHKERRQ(ierr);
193     ierr = PetscFree3(rvalues,source,recv_waits);CHKERRQ(ierr);
194     ierr = PetscFree2(nprocs,owner);CHKERRQ(ierr);
195     if (nsends) {   /* wait on sends */
196       ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
197       ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
198       ierr = PetscFree(send_status);CHKERRQ(ierr);
199     }
200     ierr = PetscFree3(svalues,send_waits,starts);CHKERRQ(ierr);
201     ierr = PetscLayoutDestroy(map);CHKERRQ(ierr);
202     ierr = PetscLayoutDestroy(nmap);CHKERRQ(ierr);
203 
204     ierr = VecCreateMPI(comm,slen,PETSC_DETERMINE,&red->b);CHKERRQ(ierr);
205     ierr = VecDuplicate(red->b,&red->x);CHKERRQ(ierr);
206     ierr = MatGetVecs(pc->pmat,&tvec,PETSC_NULL);CHKERRQ(ierr);
207     ierr = VecScatterCreate(tvec,red->is,red->b,PETSC_NULL,&red->scatter);CHKERRQ(ierr);
208     ierr = VecDestroy(tvec);CHKERRQ(ierr);
209     ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_INITIAL_MATRIX,&tmat);CHKERRQ(ierr);
210     ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
211     ierr = MatDestroy(tmat);CHKERRQ(ierr);
212   }
213 
214   /* get diagonal portion of matrix */
215   ierr = PetscMalloc(red->dcnt*sizeof(PetscScalar),&red->diag);CHKERRQ(ierr);
216   ierr = MatGetVecs(pc->pmat,&diag,PETSC_NULL);CHKERRQ(ierr);
217   ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr);
218   ierr = VecGetArrayRead(diag,&d);CHKERRQ(ierr);
219   for (i=0; i<red->dcnt; i++) {
220     red->diag[i] = 1.0/d[red->drows[i]];
221   }
222   ierr = VecRestoreArrayRead(diag,&d);CHKERRQ(ierr);
223   ierr = VecDestroy(diag);CHKERRQ(ierr);
224   ierr = KSPSetUp(red->ksp);CHKERRQ(ierr);
225   PetscFunctionReturn(0);
226 }
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "PCApply_Redistribute"
230 static PetscErrorCode PCApply_Redistribute(PC pc,Vec b,Vec x)
231 {
232   PC_Redistribute   *red = (PC_Redistribute*)pc->data;
233   PetscErrorCode    ierr;
234   PetscInt          dcnt = red->dcnt,i;
235   const PetscInt    *drows = red->drows;
236   PetscScalar       *xwork;
237   const PetscScalar *bwork,*diag = red->diag;
238 
239   PetscFunctionBegin;
240   if (!red->work) {
241     ierr = VecDuplicate(b,&red->work);CHKERRQ(ierr);
242   }
243   /* compute the rows of solution that have diagonal entries only */
244   ierr = VecSet(x,0.0);CHKERRQ(ierr);         /* x = diag(A)^{-1} b */
245   ierr = VecGetArray(x,&xwork);CHKERRQ(ierr);
246   ierr = VecGetArrayRead(b,&bwork);CHKERRQ(ierr);
247   for (i=0; i<dcnt; i++) {
248     xwork[drows[i]] = diag[i]*bwork[drows[i]];
249   }
250   ierr = PetscLogFlops(dcnt);CHKERRQ(ierr);
251   ierr = VecRestoreArray(red->work,&xwork);CHKERRQ(ierr);
252   ierr = VecRestoreArrayRead(b,&bwork);CHKERRQ(ierr);
253   /* update the right hand side for the reduced system with diagonal rows (and corresponding columns) removed */
254   ierr = MatMult(pc->pmat,x,red->work);CHKERRQ(ierr);
255   ierr = VecAYPX(red->work,-1.0,b);CHKERRQ(ierr);   /* red->work = b - A x */
256 
257   ierr = VecScatterBegin(red->scatter,red->work,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
258   ierr = VecScatterEnd(red->scatter,red->work,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
259   ierr = KSPSolve(red->ksp,red->b,red->x);CHKERRQ(ierr);
260   ierr = VecScatterBegin(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
261   ierr = VecScatterEnd(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "PCDestroy_Redistribute"
267 static PetscErrorCode PCDestroy_Redistribute(PC pc)
268 {
269   PC_Redistribute *red = (PC_Redistribute*)pc->data;
270   PetscErrorCode  ierr;
271 
272   PetscFunctionBegin;
273   if (red->scatter)  {ierr = VecScatterDestroy(red->scatter);CHKERRQ(ierr);}
274   if (red->is)       {ierr = ISDestroy(red->is);CHKERRQ(ierr);}
275   if (red->b)        {ierr = VecDestroy(red->b);CHKERRQ(ierr);}
276   if (red->x)        {ierr = VecDestroy(red->x);CHKERRQ(ierr);}
277   if (red->ksp)      {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);}
278   if (red->work)     {ierr = VecDestroy(red->work);CHKERRQ(ierr);}
279   ierr = PetscFree(red->drows);
280   ierr = PetscFree(red->diag);
281   ierr = PetscFree(red);CHKERRQ(ierr);
282   PetscFunctionReturn(0);
283 }
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "PCSetFromOptions_Redistribute"
287 static PetscErrorCode PCSetFromOptions_Redistribute(PC pc)
288 {
289   PetscErrorCode  ierr;
290   PC_Redistribute *red = (PC_Redistribute*)pc->data;
291 
292   PetscFunctionBegin;
293   ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 
297 /* -------------------------------------------------------------------------------------*/
298 /*MC
299      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
300 
301      Options for the redistribute preconditioners can be set with -redistribute_ksp_xxx <values> and -redistribute_pc_xxx <values>
302 
303      Notes:  Usually run this with -ksp_type preonly
304 
305      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
306      -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc to take advantage of the symmetry.
307 
308    Level: intermediate
309 
310 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
311 M*/
312 
313 EXTERN_C_BEGIN
314 #undef __FUNCT__
315 #define __FUNCT__ "PCCreate_Redistribute"
316 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redistribute(PC pc)
317 {
318   PetscErrorCode  ierr;
319   PC_Redistribute *red;
320   const char      *prefix;
321 
322   PetscFunctionBegin;
323   ierr = PetscNewLog(pc,PC_Redistribute,&red);CHKERRQ(ierr);
324   pc->data            = (void*)red;
325 
326   pc->ops->apply           = PCApply_Redistribute;
327   pc->ops->applytranspose  = 0;
328   pc->ops->setup           = PCSetUp_Redistribute;
329   pc->ops->destroy         = PCDestroy_Redistribute;
330   pc->ops->setfromoptions  = PCSetFromOptions_Redistribute;
331   pc->ops->view            = PCView_Redistribute;
332 
333   ierr = KSPCreate(((PetscObject)pc)->comm,&red->ksp);CHKERRQ(ierr);
334   ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
335   ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr);
336   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
337   ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
338   ierr = KSPAppendOptionsPrefix(red->ksp,"redistribute_");CHKERRQ(ierr);
339   PetscFunctionReturn(0);
340 }
341 EXTERN_C_END
342