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