xref: /petsc/src/ksp/pc/impls/redistribute/redistribute.c (revision 9b48462b03c6cfab2f47579726436d8616632354)
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   Mat          mat;
13   VecScatter   scatter;
14   IS           is;
15 } PC_Redistribute;
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "PCView_Redistribute"
19 static PetscErrorCode PCView_Redistribute(PC pc,PetscViewer viewer)
20 {
21   PC_Redistribute *red = (PC_Redistribute*)pc->data;
22   PetscErrorCode  ierr;
23   PetscTruth      iascii,isstring;
24 
25   PetscFunctionBegin;
26   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
27   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
28   if (iascii) {
29     ierr = PetscViewerASCIIPrintf(viewer,"  Redistribute preconditioner: \n");CHKERRQ(ierr);
30     ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr);
31   } else if (isstring) {
32     ierr = PetscViewerStringSPrintf(viewer," Redistribute preconditioner");CHKERRQ(ierr);
33     ierr = KSPView(red->ksp,viewer);CHKERRQ(ierr);
34   } else {
35     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redistribute",((PetscObject)viewer)->type_name);
36   }
37   PetscFunctionReturn(0);
38 }
39 
40 #include "private/matimpl.h"        /*I "petscmat.h" I*/
41 #undef __FUNCT__
42 #define __FUNCT__ "PCSetUp_Redistribute"
43 static PetscErrorCode PCSetUp_Redistribute(PC pc)
44 {
45   PC_Redistribute   *red = (PC_Redistribute*)pc->data;
46   PetscErrorCode    ierr;
47   MPI_Comm          comm;
48   PetscInt          rstart,rend,i,nz,cnt,*rows,ncnt;
49   PetscMap          *map,*nmap;
50   PetscMPIInt       size,rank,imdex,tag,n;
51   PetscInt          *source = PETSC_NULL;
52   PetscMPIInt       *nprocs = PETSC_NULL,nrecvs;
53   PetscInt          j,nsends;
54   PetscInt          *owner = PETSC_NULL,*starts = PETSC_NULL,count,slen;
55   PetscInt          *rvalues,*svalues,recvtotal;
56   PetscMPIInt       *onodes1,*olengths1;
57   MPI_Request       *send_waits = PETSC_NULL,*recv_waits = PETSC_NULL;
58   MPI_Status        recv_status,*send_status;
59   Vec               tvec;
60   Mat               tmat;
61 
62   PetscFunctionBegin;
63   if (pc->setupcalled) {
64     ierr = KSPGetOperators(red->ksp,PETSC_NULL,&tmat,PETSC_NULL);CHKERRQ(ierr);
65     ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_REUSE_MATRIX,&tmat);CHKERRQ(ierr);
66   } else {
67     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
68     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
69     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
70     ierr = PetscObjectGetNewTag((PetscObject)pc,&tag);CHKERRQ(ierr);
71 
72     /* count non-diagonal rows on process */
73     ierr = MatGetOwnershipRange(pc->mat,&rstart,&rend);CHKERRQ(ierr);
74     cnt  = 0;
75     for (i=rstart; i<rend; i++) {
76       ierr = MatGetRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
77       if (nz > 1) cnt++;
78       ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
79     }
80     ierr = PetscMalloc(cnt*sizeof(PetscInt),&rows);CHKERRQ(ierr);
81 
82     /* list non-diagonal rows on process */
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) rows[cnt++] = i;
87       ierr = MatRestoreRow(pc->mat,i,&nz,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
88     }
89 
90     /* create PetscMap for non-diagonal rows on each process */
91     ierr = PetscMalloc(sizeof(PetscMap),&map);CHKERRQ(ierr);
92     ierr = PetscMapInitialize(comm,map);CHKERRQ(ierr);
93     ierr = PetscMapSetLocalSize(map,cnt);CHKERRQ(ierr);
94     ierr = PetscMapSetBlockSize(map,1);CHKERRQ(ierr);
95     ierr = PetscMapSetUp(map);CHKERRQ(ierr);
96     rstart = map->rstart;
97     rend   = map->rend;
98 
99     /* create PetscMap for load-balanced non-diagonal rows on each process */
100     ierr = PetscMalloc(sizeof(PetscMap),&nmap);CHKERRQ(ierr);
101     ierr = PetscMapInitialize(comm,nmap);CHKERRQ(ierr);
102     ierr = MPI_Allreduce(&cnt,&ncnt,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
103     ierr = PetscMapSetSize(nmap,ncnt);CHKERRQ(ierr);
104     ierr = PetscMapSetBlockSize(nmap,1);CHKERRQ(ierr);
105     ierr = PetscMapSetUp(nmap);CHKERRQ(ierr);
106 
107     /*
108 	this code is taken from VecScatterCreate_PtoS()
109 	Determines what rows need to be moved where to
110 	load balance the non-diagonal rows
111     */
112     /*  count number of contributors to each processor */
113     ierr = PetscMalloc2(size,PetscMPIInt,&nprocs,cnt,PetscInt,&owner);CHKERRQ(ierr);
114     ierr = PetscMemzero(nprocs,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
115     j      = 0;
116     nsends = 0;
117     for (i=rstart; i<rend; i++) {
118       if (i < nmap->range[j]) j = 0;
119       for (; j<size; j++) {
120 	if (i < nmap->range[j+1]) {
121 	  if (!nprocs[j]++) nsends++;
122 	  owner[i-rstart] = j;
123 	  break;
124 	}
125       }
126     }
127     /* inform other processors of number of messages and max length*/
128     ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,nprocs,&nrecvs);CHKERRQ(ierr);
129     ierr = PetscGatherMessageLengths(comm,nsends,nrecvs,nprocs,&onodes1,&olengths1);CHKERRQ(ierr);
130     ierr = PetscSortMPIIntWithArray(nrecvs,onodes1,olengths1);CHKERRQ(ierr);
131     recvtotal = 0; for (i=0; i<nrecvs; i++) recvtotal += olengths1[i];
132 
133     /* post receives:  rvalues - rows I will own; count - nu */
134     ierr = PetscMalloc3(recvtotal,PetscInt,&rvalues,nrecvs,PetscInt,&source,nrecvs,MPI_Request,&recv_waits);CHKERRQ(ierr);
135     count  = 0;
136     for (i=0; i<nrecvs; i++) {
137       ierr  = MPI_Irecv((rvalues+count),olengths1[i],MPIU_INT,onodes1[i],tag,comm,recv_waits+i);CHKERRQ(ierr);
138       count += olengths1[i];
139     }
140 
141     /* do sends:
142        1) starts[i] gives the starting index in svalues for stuff going to
143        the ith processor
144     */
145     ierr = PetscMalloc3(cnt,PetscInt,&svalues,nsends,MPI_Request,&send_waits,size,PetscInt,&starts);CHKERRQ(ierr);
146     starts[0]  = 0;
147     for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
148     for (i=0; i<cnt; i++) {
149       svalues[starts[owner[i]]++] = rows[i];
150     }
151     ierr = PetscFree(rows);CHKERRQ(ierr);
152     starts[0] = 0;
153     for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[i-1];}
154     count = 0;
155     for (i=0; i<size; i++) {
156       if (nprocs[i]) {
157 	ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
158       }
159     }
160 
161     /*  wait on receives */
162     count  = nrecvs;
163     slen   = 0;
164     while (count) {
165       ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
166       /* unpack receives into our local space */
167       ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
168       slen += n;
169       count--;
170     }
171     if (slen != recvtotal) SETERRQ2(PETSC_ERR_PLIB,"Total message lengths %D not expected %D",slen,recvtotal);
172 
173     ierr = ISCreateGeneral(comm,slen,rvalues,&red->is);CHKERRQ(ierr);
174 
175     /* free up all work space */
176     ierr = PetscFree(olengths1);CHKERRQ(ierr);
177     ierr = PetscFree(onodes1);CHKERRQ(ierr);
178     ierr = PetscFree3(rvalues,source,recv_waits);CHKERRQ(ierr);
179     ierr = PetscFree2(nprocs,owner);CHKERRQ(ierr);
180     if (nsends) {   /* wait on sends */
181       ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
182       ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
183       ierr = PetscFree(send_status);CHKERRQ(ierr);
184     }
185     ierr = PetscFree3(svalues,send_waits,starts);CHKERRQ(ierr);
186     ierr = PetscFree(map->range);CHKERRQ(ierr);
187     ierr = PetscFree(map);CHKERRQ(ierr);
188     ierr = PetscFree(nmap->range);CHKERRQ(ierr);
189     ierr = PetscFree(nmap);CHKERRQ(ierr);
190 
191     ierr = VecCreateMPI(comm,slen,PETSC_DETERMINE,&red->b);CHKERRQ(ierr);
192     ierr = VecDuplicate(red->b,&red->x);CHKERRQ(ierr);
193     ierr = MatGetVecs(pc->pmat,&tvec,PETSC_NULL);CHKERRQ(ierr);
194     ierr = VecScatterCreate(tvec,red->is,red->b,PETSC_NULL,&red->scatter);CHKERRQ(ierr);
195     ierr = VecDestroy(tvec);CHKERRQ(ierr);
196     ierr = MatGetSubMatrix(pc->pmat,red->is,red->is,MAT_INITIAL_MATRIX,&tmat);CHKERRQ(ierr);
197   }
198   ierr = KSPSetOperators(red->ksp,tmat,tmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
199   ierr = MatDestroy(tmat);CHKERRQ(ierr);
200   ierr = KSPSetUp(red->ksp);CHKERRQ(ierr);
201   PetscFunctionReturn(0);
202 }
203 
204 #undef __FUNCT__
205 #define __FUNCT__ "PCApply_Redistribute"
206 static PetscErrorCode PCApply_Redistribute(PC pc,Vec b,Vec x)
207 {
208   PC_Redistribute   *red = (PC_Redistribute*)pc->data;
209   PetscErrorCode    ierr;
210 
211   PetscFunctionBegin;
212   ierr = VecScatterBegin(red->scatter,b,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
213   ierr = VecScatterEnd(red->scatter,b,red->b,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
214   ierr = KSPSolve(red->ksp,red->b,red->x);CHKERRQ(ierr);
215   ierr = VecScatterBegin(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
216   ierr = VecScatterEnd(red->scatter,red->x,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
217   PetscFunctionReturn(0);
218 }
219 
220 #undef __FUNCT__
221 #define __FUNCT__ "PCDestroy_Redistribute"
222 static PetscErrorCode PCDestroy_Redistribute(PC pc)
223 {
224   PC_Redistribute *red = (PC_Redistribute*)pc->data;
225   PetscErrorCode  ierr;
226 
227   PetscFunctionBegin;
228   if (red->scatter)  {ierr = VecScatterDestroy(red->scatter);CHKERRQ(ierr);}
229   if (red->is)       {ierr = ISDestroy(red->is);CHKERRQ(ierr);}
230   if (red->b)        {ierr = VecDestroy(red->b);CHKERRQ(ierr);}
231   if (red->x)        {ierr = VecDestroy(red->x);CHKERRQ(ierr);}
232   if (red->mat)      {ierr = MatDestroy(red->mat);CHKERRQ(ierr);}
233   if (red->ksp)      {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);}
234   ierr = PetscFree(red);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "PCSetFromOptions_Redistribute"
240 static PetscErrorCode PCSetFromOptions_Redistribute(PC pc)
241 {
242   PetscErrorCode  ierr;
243   PC_Redistribute *red = (PC_Redistribute*)pc->data;
244 
245   PetscFunctionBegin;
246   ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr);
247   PetscFunctionReturn(0);
248 }
249 
250 /* -------------------------------------------------------------------------------------*/
251 /*MC
252      PCREDISTRIBUTE - Redistributes a matrix for load balancing and then applys a KSP to that new matrix
253 
254      Options for the redistribute preconditioners can be set with -redistribute_ksp_xxx <values> and -redistribute_pc_xxx <values>
255 
256    Level: intermediate
257 
258 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
259 M*/
260 
261 EXTERN_C_BEGIN
262 #undef __FUNCT__
263 #define __FUNCT__ "PCCreate_Redistribute"
264 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redistribute(PC pc)
265 {
266   PetscErrorCode  ierr;
267   PC_Redistribute *red;
268   const char      *prefix;
269 
270   PetscFunctionBegin;
271   ierr = PetscNewLog(pc,PC_Redistribute,&red);CHKERRQ(ierr);
272   pc->data            = (void*)red;
273 
274   pc->ops->apply           = PCApply_Redistribute;
275   pc->ops->applytranspose  = 0;
276   pc->ops->setup           = PCSetUp_Redistribute;
277   pc->ops->destroy         = PCDestroy_Redistribute;
278   pc->ops->setfromoptions  = PCSetFromOptions_Redistribute;
279   pc->ops->view            = PCView_Redistribute;
280 
281   ierr = KSPCreate(((PetscObject)pc)->comm,&red->ksp);CHKERRQ(ierr);
282   ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
283   ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr);
284   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
285   ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
286   ierr = KSPAppendOptionsPrefix(red->ksp,"redistribute_");CHKERRQ(ierr);
287   PetscFunctionReturn(0);
288 }
289 EXTERN_C_END
290