xref: /petsc/src/ksp/pc/impls/gamg/geo.c (revision fee518eb64eefbf9bda05b76e476d2dcf413f412)
1 /*
2  GAMG geometric-algebric multiogrid PC - Mark Adams 2011
3  */
4 
5 #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
6 #include <petsc-private/kspimpl.h>
7 
8 #if defined(PETSC_HAVE_TRIANGLE)
9 #define REAL PetscReal
10 #include <triangle.h>
11 #endif
12 
13 #include <petscblaslapack.h>
14 
15 /* Private context for the GAMG preconditioner */
16 typedef struct {
17   PetscInt lid;            /* local vertex index */
18   PetscInt degree;         /* vertex degree */
19 } GAMGNode;
20 
21 int petsc_geo_mg_compare(const void *a, const void *b)
22 {
23   return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree);
24 }
25 
26 /* -------------------------------------------------------------------------- */
27 /*
28    PCSetCoordinates_GEO
29 
30    Input Parameter:
31    .  pc - the preconditioner context
32 */
33 EXTERN_C_BEGIN
34 #undef __FUNCT__
35 #define __FUNCT__ "PCSetCoordinates_GEO"
36 PetscErrorCode PCSetCoordinates_GEO(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
37 {
38   PC_MG          *mg      = (PC_MG*)pc->data;
39   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
40   PetscErrorCode ierr;
41   PetscInt       arrsz,bs,my0,kk,ii,nloc,Iend;
42   Mat            Amat = pc->pmat;
43 
44   PetscFunctionBegin;
45   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
46   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
47 
48   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
49   nloc = (Iend-my0)/bs;
50 
51   if (nloc!=a_nloc) SETERRQ2(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Stokes not supported nloc = %d %d.",a_nloc,nloc);
52   if ((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
53 
54   pc_gamg->data_cell_rows = 1;
55   if (coords==0 && nloc > 0) SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'.");
56   pc_gamg->data_cell_cols = ndm; /* coordinates */
57 
58   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
59 
60   /* create data - syntactic sugar that should be refactored at some point */
61   if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) {
62     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
63     ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr);
64   }
65   for (kk=0; kk<arrsz; kk++) pc_gamg->data[kk] = -999.;
66   pc_gamg->data[arrsz] = -99.;
67   /* copy data in - column oriented */
68   for (kk = 0; kk < nloc; kk++) {
69     for (ii = 0; ii < ndm; ii++) {
70       pc_gamg->data[ii*nloc + kk] =  coords[kk*ndm + ii];
71     }
72   }
73   if (pc_gamg->data[arrsz] != -99.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data[arrsz %D] %g != -99.",arrsz,pc_gamg->data[arrsz]);
74   pc_gamg->data_sz = arrsz;
75   PetscFunctionReturn(0);
76 }
77 EXTERN_C_END
78 
79 /* -------------------------------------------------------------------------- */
80 /*
81    PCSetData_GEO
82 
83   Input Parameter:
84    . pc -
85 */
86 #undef __FUNCT__
87 #define __FUNCT__ "PCSetData_GEO"
88 PetscErrorCode PCSetData_GEO(PC pc, Mat m)
89 {
90   PetscFunctionBegin;
91   SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"GEO MG needs coordinates");
92 }
93 
94 /* -------------------------------------------------------------------------- */
95 /*
96    PCSetFromOptions_GEO
97 
98   Input Parameter:
99    . pc -
100 */
101 #undef __FUNCT__
102 #define __FUNCT__ "PCSetFromOptions_GEO"
103 PetscErrorCode PCSetFromOptions_GEO(PC pc)
104 {
105   PetscErrorCode ierr;
106   PC_MG          *mg      = (PC_MG*)pc->data;
107   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
108 
109   PetscFunctionBegin;
110   ierr = PetscOptionsHead("GAMG-GEO options");CHKERRQ(ierr);
111   {
112     /* -pc_gamg_sa_nsmooths */
113     /* pc_gamg_sa->smooths = 0; */
114     /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */
115     /*                        "smoothing steps for smoothed aggregation, usually 1 (0)", */
116     /*                        "PCGAMGSetNSmooths_AGG", */
117     /*                        pc_gamg_sa->smooths, */
118     /*                        &pc_gamg_sa->smooths, */
119     /*                        &flag);  */
120     /* CHKERRQ(ierr); */
121   }
122   ierr = PetscOptionsTail();CHKERRQ(ierr);
123 
124   /* call base class */
125   ierr = PCSetFromOptions_GAMG(pc);CHKERRQ(ierr);
126 
127   if (pc_gamg->verbose) {
128     MPI_Comm wcomm = ((PetscObject)pc)->comm;
129     PetscPrintf(wcomm,"[%d]%s done\n",0,__FUNCT__);
130   }
131   PetscFunctionReturn(0);
132 }
133 
134 /* -------------------------------------------------------------------------- */
135 /*
136  triangulateAndFormProl
137 
138    Input Parameter:
139    . selected_2 - list of selected local ID, includes selected ghosts
140    . data_stride -
141    . coords[2*data_stride] - column vector of local coordinates w/ ghosts
142    . nselected_1 - selected IDs that go with base (1) graph
143    . clid_lid_1[nselected_1] - lids of selected (c) nodes   ???????????
144    . agg_lists_1 - list of aggregates
145    . crsGID[selected.size()] - global index for prolongation operator
146    . bs - block size
147   Output Parameter:
148    . a_Prol - prolongation operator
149    . a_worst_best - measure of worst missed fine vertex, 0 is no misses
150 */
151 #undef __FUNCT__
152 #define __FUNCT__ "triangulateAndFormProl"
153 static PetscErrorCode triangulateAndFormProl(IS  selected_2, /* list of selected local ID, includes selected ghosts */
154                                              const PetscInt data_stride,
155                                              const PetscReal coords[], /* column vector of local coordinates w/ ghosts */
156                                              const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */
157                                              const PetscInt clid_lid_1[],
158                                              const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */
159                                              const PetscInt crsGID[],
160                                              const PetscInt bs,
161                                              Mat a_Prol, /* prolongation operator (output) */
162                                              PetscReal *a_worst_best) /* measure of worst missed fine vertex, 0 is no misses */
163 {
164 #if defined(PETSC_HAVE_TRIANGLE)
165   PetscErrorCode       ierr;
166   PetscInt             jj,tid,tt,idx,nselected_2;
167   struct triangulateio in,mid;
168   const PetscInt       *selected_idx_2;
169   PetscMPIInt          rank,size;
170   PetscInt             Istart,Iend,nFineLoc,myFine0;
171   int                  kk,nPlotPts,sid;
172   MPI_Comm             wcomm = ((PetscObject)a_Prol)->comm;
173   PetscReal            tm;
174 
175   PetscFunctionBegin;
176   ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
177   ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
178   ierr = ISGetSize(selected_2, &nselected_2);CHKERRQ(ierr);
179   if (nselected_2 == 1 || nselected_2 == 2) { /* 0 happens on idle processors */
180     *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */
181   } else *a_worst_best = 0.0;
182   ierr = MPI_Allreduce(a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr);
183   if (tm > 0.0) {
184     *a_worst_best = 100.0;
185     PetscFunctionReturn(0);
186   }
187   ierr     = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
188   nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs;
189   nPlotPts = nFineLoc; /* locals */
190   /* traingle */
191   /* Define input points - in*/
192   in.numberofpoints          = nselected_2;
193   in.numberofpointattributes = 0;
194   /* get nselected points */
195   ierr = PetscMalloc(2*(nselected_2)*sizeof(REAL), &in.pointlist);CHKERRQ(ierr);
196   ierr = ISGetIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
197 
198   for (kk=0,sid=0; kk<nselected_2; kk++,sid += 2) {
199     PetscInt lid = selected_idx_2[kk];
200     in.pointlist[sid]   = coords[lid];
201     in.pointlist[sid+1] = coords[data_stride + lid];
202     if (lid>=nFineLoc) nPlotPts++;
203   }
204   if (sid != 2*nselected_2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != 2*nselected_2 %D",sid,nselected_2);
205 
206   in.numberofsegments      = 0;
207   in.numberofedges         = 0;
208   in.numberofholes         = 0;
209   in.numberofregions       = 0;
210   in.trianglelist          = 0;
211   in.segmentmarkerlist     = 0;
212   in.pointattributelist    = 0;
213   in.pointmarkerlist       = 0;
214   in.triangleattributelist = 0;
215   in.trianglearealist      = 0;
216   in.segmentlist           = 0;
217   in.holelist              = 0;
218   in.regionlist            = 0;
219   in.edgelist              = 0;
220   in.edgemarkerlist        = 0;
221   in.normlist              = 0;
222 
223   /* triangulate */
224   mid.pointlist = 0;            /* Not needed if -N switch used. */
225   /* Not needed if -N switch used or number of point attributes is zero: */
226   mid.pointattributelist = 0;
227   mid.pointmarkerlist    = 0; /* Not needed if -N or -B switch used. */
228   mid.trianglelist       = 0;    /* Not needed if -E switch used. */
229   /* Not needed if -E switch used or number of triangle attributes is zero: */
230   mid.triangleattributelist = 0;
231   mid.neighborlist          = 0; /* Needed only if -n switch used. */
232   /* Needed only if segments are output (-p or -c) and -P not used: */
233   mid.segmentlist = 0;
234   /* Needed only if segments are output (-p or -c) and -P and -B not used: */
235   mid.segmentmarkerlist = 0;
236   mid.edgelist          = 0;    /* Needed only if -e switch used. */
237   mid.edgemarkerlist    = 0; /* Needed if -e used and -B not used. */
238   mid.numberoftriangles = 0;
239 
240   /* Triangulate the points.  Switches are chosen to read and write a  */
241   /*   PSLG (p), preserve the convex hull (c), number everything from  */
242   /*   zero (z), assign a regional attribute to each element (A), and  */
243   /*   produce an edge list (e), a Voronoi diagram (v), and a triangle */
244   /*   neighbor list (n).                                            */
245   if (nselected_2 != 0) { /* inactive processor */
246     char args[] = "npczQ"; /* c is needed ? */
247     triangulate(args, &in, &mid, (struct triangulateio*) NULL);
248     /* output .poly files for 'showme' */
249     if (!PETSC_TRUE) {
250       static int level = 1;
251       FILE       *file; char fname[32];
252 
253       sprintf(fname,"C%d_%d.poly",level,rank); file = fopen(fname, "w");
254       /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/
255       fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0);
256       /*Following lines: <vertex #> <x> <y> */
257       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid += 2) {
258         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
259       }
260       /*One line: <# of segments> <# of boundary markers (0 or 1)> */
261       fprintf(file, "%d  %d\n",0,0);
262       /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */
263       /* One line: <# of holes> */
264       fprintf(file, "%d\n",0);
265       /* Following lines: <hole #> <x> <y> */
266       /* Optional line: <# of regional attributes and/or area constraints> */
267       /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */
268       fclose(file);
269 
270       /* elems */
271       sprintf(fname,"C%d_%d.ele",level,rank); file = fopen(fname, "w");
272       /* First line: <# of triangles> <nodes per triangle> <# of attributes> */
273       fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0);
274       /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */
275       for (kk=0,sid=0; kk<mid.numberoftriangles; kk++,sid += 3) {
276         fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]);
277       }
278       fclose(file);
279 
280       sprintf(fname,"C%d_%d.node",level,rank); file = fopen(fname, "w");
281       /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */
282       /* fprintf(file, "%d  %d  %d  %d\n",in.numberofpoints,2,0,0); */
283       fprintf(file, "%d  %d  %d  %d\n",nPlotPts,2,0,0);
284       /*Following lines: <vertex #> <x> <y> */
285       for (kk=0,sid=0; kk<in.numberofpoints; kk++,sid+=2) {
286         fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]);
287       }
288 
289       sid /= 2;
290       for (jj=0; jj<nFineLoc; jj++) {
291         PetscBool sel = PETSC_TRUE;
292         for (kk=0; kk<nselected_2 && sel; kk++) {
293           PetscInt lid = selected_idx_2[kk];
294           if (lid == jj) sel = PETSC_FALSE;
295         }
296         if (sel) fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[data_stride + jj]);
297       }
298       fclose(file);
299       if (sid != nPlotPts) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sid %D != nPlotPts %D",sid,nPlotPts);
300       level++;
301     }
302   }
303 #if defined PETSC_GAMG_USE_LOG
304   ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
305 #endif
306   { /* form P - setup some maps */
307     PetscInt clid,mm,*nTri,*node_tri;
308 
309     ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &node_tri);CHKERRQ(ierr);
310     ierr = PetscMalloc(nselected_2*sizeof(PetscInt), &nTri);CHKERRQ(ierr);
311 
312     /* need list of triangles on node */
313     for (kk=0; kk<nselected_2; kk++) nTri[kk] = 0;
314     for (tid=0,kk=0; tid<mid.numberoftriangles; tid++) {
315       for (jj=0; jj<3; jj++) {
316         PetscInt cid = mid.trianglelist[kk++];
317         if (nTri[cid] == 0) node_tri[cid] = tid;
318         nTri[cid]++;
319       }
320     }
321 #define EPS 1.e-12
322     /* find points and set prolongation */
323     for (mm = clid = 0; mm < nFineLoc; mm++) {
324       PetscBool ise;
325       ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise);CHKERRQ(ierr);
326       if (!ise) {
327         const PetscInt lid = mm;
328         /* for (clid_iterator=0;clid_iterator<nselected_1;clid_iterator++) { */
329         PetscScalar  AA[3][3];
330         PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO;
331         PetscCDPos   pos;
332         ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
333         while (pos) {
334           PetscInt flid;
335           ierr = PetscLLNGetID(pos, &flid);CHKERRQ(ierr);
336           ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos);CHKERRQ(ierr);
337 
338           if (flid < nFineLoc) {  /* could be a ghost */
339             PetscInt       bestTID = -1; PetscReal best_alpha = 1.e10;
340             const PetscInt fgid    = flid + myFine0;
341             /* compute shape function for gid */
342             const PetscReal fcoord[3] = {coords[flid],coords[data_stride+flid],1.0};
343             PetscBool       haveit    =PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3];
344 
345             /* look for it */
346             for (tid = node_tri[clid], jj=0;
347                  jj < 5 && !haveit && tid != -1;
348                  jj++) {
349               for (tt=0; tt<3; tt++) {
350                 PetscInt cid2 = mid.trianglelist[3*tid + tt];
351                 PetscInt lid2 = selected_idx_2[cid2];
352                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
353                 clids[tt] = cid2; /* store for interp */
354               }
355 
356               for (tt=0; tt<3; tt++) alpha[tt] = (PetscScalar)fcoord[tt];
357 
358               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
359               PetscStackCall("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
360               {
361                 PetscBool have=PETSC_TRUE;  PetscReal lowest=1.e10;
362                 for (tt = 0, idx = 0; tt < 3; tt++) {
363                   if (PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS) have = PETSC_FALSE;
364                   if (PetscRealPart(alpha[tt]) < lowest) {
365                     lowest = PetscRealPart(alpha[tt]);
366                     idx    = tt;
367                   }
368                 }
369                 haveit = have;
370               }
371               tid = mid.neighborlist[3*tid + idx];
372             }
373 
374             if (!haveit) {
375               /* brute force */
376               for (tid=0; tid<mid.numberoftriangles && !haveit; tid++) {
377                 for (tt=0; tt<3; tt++) {
378                   PetscInt cid2 = mid.trianglelist[3*tid + tt];
379                   PetscInt lid2 = selected_idx_2[cid2];
380                   AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
381                   clids[tt] = cid2; /* store for interp */
382                 }
383                 for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
384                 /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
385                 PetscStackCall("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
386                 {
387                   PetscBool have=PETSC_TRUE;  PetscReal worst=0.0, v;
388                   for (tt=0; tt<3 && have; tt++) {
389                     if (PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS) have=PETSC_FALSE;
390                     if ((v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst) worst = v;
391                   }
392                   if (worst < best_alpha) {
393                     best_alpha = worst; bestTID = tid;
394                   }
395                   haveit = have;
396                 }
397               }
398             }
399             if (!haveit) {
400               if (best_alpha > *a_worst_best) *a_worst_best = best_alpha;
401               /* use best one */
402               for (tt=0; tt<3; tt++) {
403                 PetscInt cid2 = mid.trianglelist[3*bestTID + tt];
404                 PetscInt lid2 = selected_idx_2[cid2];
405                 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[data_stride + lid2]; AA[tt][2] = 1.0;
406                 clids[tt] = cid2; /* store for interp */
407               }
408               for (tt=0; tt<3; tt++) alpha[tt] = fcoord[tt];
409               /* SUBROUTINE DGESV(N, NRHS, A, LDA, IPIV, B, LDB, INFO) */
410               PetscStackCall("LAPACKgesv",LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO));
411             }
412 
413             /* put in row of P */
414             for (idx=0; idx<3; idx++) {
415               PetscScalar shp = alpha[idx];
416               if (PetscAbs(PetscRealPart(shp)) > 1.e-6) {
417                 PetscInt cgid = crsGID[clids[idx]];
418                 PetscInt jj   = cgid*bs, ii = fgid*bs; /* need to gloalize */
419                 for (tt=0; tt < bs; tt++, ii++, jj++) {
420                   ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES);CHKERRQ(ierr);
421                 }
422               }
423             }
424           }
425         } /* aggregates iterations */
426         clid++;
427       } /* a coarse agg */
428     } /* for all fine nodes */
429 
430     ierr = ISRestoreIndices(selected_2, &selected_idx_2);CHKERRQ(ierr);
431     ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
432     ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
433 
434     ierr = PetscFree(node_tri);CHKERRQ(ierr);
435     ierr = PetscFree(nTri);CHKERRQ(ierr);
436   }
437 #if defined PETSC_GAMG_USE_LOG
438   ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr);
439 #endif
440   free(mid.trianglelist);
441   free(mid.neighborlist);
442   ierr = PetscFree(in.pointlist);CHKERRQ(ierr);
443   PetscFunctionReturn(0);
444 #else
445   SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_PLIB,"configure with TRIANGLE to use geometric MG");
446 #endif
447 }
448 /* -------------------------------------------------------------------------- */
449 /*
450    getGIDsOnSquareGraph - square graph, get
451 
452    Input Parameter:
453    . nselected_1 - selected local indices (includes ghosts in input Gmat1)
454    . clid_lid_1 - [nselected_1] lids of selected nodes
455    . Gmat1 - graph that goes with 'selected_1'
456    Output Parameter:
457    . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2)
458    . a_Gmat_2 - graph that is squared of 'Gmat_1'
459    . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes
460 */
461 #undef __FUNCT__
462 #define __FUNCT__ "getGIDsOnSquareGraph"
463 static PetscErrorCode getGIDsOnSquareGraph(const PetscInt nselected_1,const PetscInt clid_lid_1[],const Mat Gmat1,IS *a_selected_2,Mat *a_Gmat_2,PetscInt **a_crsGID)
464 {
465   PetscErrorCode ierr;
466   PetscMPIInt    rank,size;
467   PetscInt       *crsGID, kk,my0,Iend,nloc;
468   MPI_Comm       wcomm = ((PetscObject)Gmat1)->comm;
469 
470   PetscFunctionBegin;
471   ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
472   ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
473   ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend);CHKERRQ(ierr); /* AIJ */
474   nloc = Iend - my0; /* this does not change */
475 
476   if (size == 1) { /* not much to do in serial */
477     ierr = PetscMalloc(nselected_1*sizeof(PetscInt), &crsGID);CHKERRQ(ierr);
478     for (kk=0; kk<nselected_1; kk++) crsGID[kk] = kk;
479     *a_Gmat_2 = 0;
480     ierr      = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2);CHKERRQ(ierr);
481   } else {
482     PetscInt    idx,num_fine_ghosts,num_crs_ghost,myCrs0;
483     Mat_MPIAIJ  *mpimat2;
484     Mat         Gmat2;
485     Vec         locState;
486     PetscScalar *cpcol_state;
487 
488     /* scan my coarse zero gid, set 'lid_state' with coarse GID */
489     kk = nselected_1;
490     MPI_Scan(&kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm);
491     myCrs0 -= nselected_1;
492 
493     if (a_Gmat_2) { /* output */
494       /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */
495       ierr      = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
496       *a_Gmat_2 = Gmat2; /* output */
497     } else Gmat2 = Gmat1;  /* use local to get crsGIDs at least */
498     /* get coarse grid GIDS for selected (locals and ghosts) */
499     mpimat2 = (Mat_MPIAIJ*)Gmat2->data;
500     ierr    = MatGetVecs(Gmat2, &locState, 0);CHKERRQ(ierr);
501     ierr    = VecSet(locState, (PetscScalar)(PetscReal)(-1));CHKERRQ(ierr); /* set with UNKNOWN state */
502     for (kk=0; kk<nselected_1; kk++) {
503       PetscInt    fgid = clid_lid_1[kk] + my0;
504       PetscScalar v    = (PetscScalar)(kk+myCrs0);
505       ierr = VecSetValues(locState, 1, &fgid, &v, INSERT_VALUES);CHKERRQ(ierr); /* set with PID */
506     }
507     ierr = VecAssemblyBegin(locState);CHKERRQ(ierr);
508     ierr = VecAssemblyEnd(locState);CHKERRQ(ierr);
509     ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
510     ierr =   VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
511     ierr = VecGetLocalSize(mpimat2->lvec, &num_fine_ghosts);CHKERRQ(ierr);
512     ierr = VecGetArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
513     for (kk=0,num_crs_ghost=0; kk<num_fine_ghosts; kk++) {
514       if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) num_crs_ghost++;
515     }
516     ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID);CHKERRQ(ierr); /* output */
517     {
518       PetscInt *selected_set;
519       ierr = PetscMalloc((nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set);CHKERRQ(ierr);
520       /* do ghost of 'crsGID' */
521       for (kk=0,idx=nselected_1; kk<num_fine_ghosts; kk++) {
522         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
523           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
524           selected_set[idx] = nloc + kk;
525           crsGID[idx++]     = cgid;
526         }
527       }
528       if (idx != (nselected_1+num_crs_ghost)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != (nselected_1 %D + num_crs_ghost %D)",idx,nselected_1,num_crs_ghost);
529       ierr = VecRestoreArray(mpimat2->lvec, &cpcol_state);CHKERRQ(ierr);
530       /* do locals in 'crsGID' */
531       ierr = VecGetArray(locState, &cpcol_state);CHKERRQ(ierr);
532       for (kk=0,idx=0; kk<nloc; kk++) {
533         if ((PetscInt)PetscRealPart(cpcol_state[kk]) != -1) {
534           PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]);
535           selected_set[idx] = kk;
536           crsGID[idx++]     = cgid;
537         }
538       }
539       if (idx != nselected_1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"idx %D != nselected_1 %D",idx,nselected_1);
540       ierr = VecRestoreArray(locState, &cpcol_state);CHKERRQ(ierr);
541 
542       if (a_selected_2 != 0) { /* output */
543         ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2);CHKERRQ(ierr);
544       } else {
545         ierr = PetscFree(selected_set);CHKERRQ(ierr);
546       }
547     }
548     ierr = VecDestroy(&locState);CHKERRQ(ierr);
549   }
550   *a_crsGID = crsGID; /* output */
551   PetscFunctionReturn(0);
552 }
553 
554 /* -------------------------------------------------------------------------- */
555 /*
556    PCGAMGgraph_GEO
557 
558   Input Parameter:
559    . pc - this
560    . Amat - matrix on this fine level
561   Output Parameter:
562    . a_Gmat
563 */
564 #undef __FUNCT__
565 #define __FUNCT__ "PCGAMGgraph_GEO"
566 PetscErrorCode PCGAMGgraph_GEO(PC pc,const Mat Amat,Mat *a_Gmat)
567 {
568   PetscErrorCode  ierr;
569   PC_MG           *mg      = (PC_MG*)pc->data;
570   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
571   const PetscInt  verbose  = pc_gamg->verbose;
572   const PetscReal vfilter  = pc_gamg->threshold;
573   PetscMPIInt     rank,size;
574   MPI_Comm        wcomm = ((PetscObject)Amat)->comm;
575   Mat             Gmat;
576   PetscBool       set,flg,symm;
577 
578   PetscFunctionBegin;
579 #if defined PETSC_USE_LOG
580   ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
581 #endif
582   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
583   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
584 
585   ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr);
586   symm = (PetscBool)!(set && flg);
587 
588   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
589   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr);
590 
591   *a_Gmat = Gmat;
592 #if defined PETSC_USE_LOG
593   ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr);
594 #endif
595   PetscFunctionReturn(0);
596 }
597 
598 /* -------------------------------------------------------------------------- */
599 /*
600    PCGAMGcoarsen_GEO
601 
602   Input Parameter:
603    . a_pc - this
604    . a_Gmat - graph
605   Output Parameter:
606    . a_llist_parent - linked list from selected indices for data locality only
607 */
608 #undef __FUNCT__
609 #define __FUNCT__ "PCGAMGcoarsen_GEO"
610 PetscErrorCode PCGAMGcoarsen_GEO(PC a_pc,Mat *a_Gmat,PetscCoarsenData **a_llist_parent)
611 {
612   PetscErrorCode ierr;
613   PC_MG          *mg      = (PC_MG*)a_pc->data;
614   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
615   PetscInt       Istart,Iend,nloc,kk,Ii,ncols;
616   PetscMPIInt    rank,size;
617   IS             perm;
618   GAMGNode       *gnodes;
619   PetscInt       *permute;
620   Mat            Gmat  = *a_Gmat;
621   MPI_Comm       wcomm = ((PetscObject)Gmat)->comm;
622   MatCoarsen     crs;
623 
624   PetscFunctionBegin;
625 #if defined PETSC_USE_LOG
626   ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
627 #endif
628   ierr = MPI_Comm_rank(wcomm, &rank);CHKERRQ(ierr);
629   ierr = MPI_Comm_size(wcomm, &size);CHKERRQ(ierr);
630   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
631   nloc = (Iend-Istart);
632 
633   /* create random permutation with sort for geo-mg */
634   ierr = PetscMalloc(nloc*sizeof(GAMGNode), &gnodes);CHKERRQ(ierr);
635   ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr);
636 
637   for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */
638     ierr = MatGetRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
639     {
640       PetscInt lid = Ii - Istart;
641       gnodes[lid].lid    = lid;
642       gnodes[lid].degree = ncols;
643     }
644     ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0);CHKERRQ(ierr);
645   }
646   /* randomize */
647   srand(1); /* make deterministic */
648   if (PETSC_TRUE) {
649     PetscBool *bIndexSet;
650     ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr);
651     for (Ii = 0; Ii < nloc; Ii++) bIndexSet[Ii] = PETSC_FALSE;
652     for (Ii = 0; Ii < nloc; Ii++) {
653       PetscInt iSwapIndex = rand()%nloc;
654       if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
655         GAMGNode iTemp = gnodes[iSwapIndex];
656         gnodes[iSwapIndex]    = gnodes[Ii];
657         gnodes[Ii]            = iTemp;
658         bIndexSet[Ii]         = PETSC_TRUE;
659         bIndexSet[iSwapIndex] = PETSC_TRUE;
660       }
661     }
662     ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
663   }
664   /* only sort locals */
665   qsort(gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare);
666   /* create IS of permutation */
667   for (kk=0; kk<nloc; kk++) permute[kk] = gnodes[kk].lid; /* locals only */
668   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
669 
670   ierr = PetscFree(gnodes);CHKERRQ(ierr);
671 
672   /* get MIS aggs */
673 
674   ierr = MatCoarsenCreate(wcomm, &crs);CHKERRQ(ierr);
675   ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr);
676   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
677   ierr = MatCoarsenSetAdjacency(crs, Gmat);CHKERRQ(ierr);
678   ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr);
679   ierr = MatCoarsenSetStrictAggs(crs, PETSC_FALSE);CHKERRQ(ierr);
680   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
681   ierr = MatCoarsenGetData(crs, a_llist_parent);CHKERRQ(ierr);
682   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
683 
684   ierr = ISDestroy(&perm);CHKERRQ(ierr);
685 #if defined PETSC_USE_LOG
686   ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr);
687 #endif
688   PetscFunctionReturn(0);
689 }
690 
691 /* -------------------------------------------------------------------------- */
692 /*
693  PCGAMGProlongator_GEO
694 
695  Input Parameter:
696  . pc - this
697  . Amat - matrix on this fine level
698  . Graph - used to get ghost data for nodes in
699  . selected_1 - [nselected]
700  . agg_lists - [nselected]
701  Output Parameter:
702  . a_P_out - prolongation operator to the next level
703  */
704 #undef __FUNCT__
705 #define __FUNCT__ "PCGAMGProlongator_GEO"
706 PetscErrorCode PCGAMGProlongator_GEO(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
707 {
708   PC_MG          *mg      = (PC_MG*)pc->data;
709   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
710   const PetscInt verbose  = pc_gamg->verbose;
711   const PetscInt dim      = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols;
712   PetscErrorCode ierr;
713   PetscInt       Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid;
714   Mat            Prol;
715   PetscMPIInt    rank, size;
716   MPI_Comm       wcomm = ((PetscObject)Amat)->comm;
717   IS             selected_2,selected_1;
718   const PetscInt *selected_idx;
719 
720   PetscFunctionBegin;
721 #if defined PETSC_USE_LOG
722   ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
723 #endif
724   ierr = MPI_Comm_rank(wcomm,&rank);CHKERRQ(ierr);
725   ierr = MPI_Comm_size(wcomm,&size);CHKERRQ(ierr);
726   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
727   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
728   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
729   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) % bs %D",Iend,Istart,bs);
730 
731   /* get 'nLocalSelected' */
732   ierr = PetscCDGetMIS(agg_lists, &selected_1);CHKERRQ(ierr);
733   ierr = ISGetSize(selected_1, &jj);CHKERRQ(ierr);
734   ierr = PetscMalloc(jj*sizeof(PetscInt), &clid_flid);CHKERRQ(ierr);
735   ierr = ISGetIndices(selected_1, &selected_idx);CHKERRQ(ierr);
736   for (kk=0,nLocalSelected=0; kk<jj; kk++) {
737     PetscInt lid = selected_idx[kk];
738     if (lid<nloc) {
739       ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
740       if (ncols>1) clid_flid[nLocalSelected++] = lid; /* fiter out singletons */
741       ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0);CHKERRQ(ierr);
742     }
743   }
744   ierr = ISRestoreIndices(selected_1, &selected_idx);CHKERRQ(ierr);
745   ierr = ISDestroy(&selected_1);CHKERRQ(ierr); /* this is selected_1 in serial */
746 
747   /* create prolongator, create P matrix */
748   ierr = MatCreate(wcomm, &Prol);CHKERRQ(ierr);
749   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
750   ierr = MatSetBlockSizes(Prol, bs, bs);CHKERRQ(ierr);
751   ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr);
752   ierr = MatSeqAIJSetPreallocation(Prol,3*data_cols,NULL);CHKERRQ(ierr);
753   ierr = MatMPIAIJSetPreallocation(Prol,3*data_cols,NULL,3*data_cols,NULL);CHKERRQ(ierr);
754   /* ierr = MatCreateAIJ(wcomm,  */
755   /*                      nloc*bs, nLocalSelected*bs, */
756   /*                      PETSC_DETERMINE, PETSC_DETERMINE, */
757   /*                      3*data_cols, NULL,  */
758   /*                      3*data_cols, NULL, */
759   /*                      &Prol); */
760   /* CHKERRQ(ierr); */
761 
762   /* can get all points "removed" - but not on geomg */
763   ierr =  MatGetSize(Prol, &kk, &jj);CHKERRQ(ierr);
764   if (jj==0) {
765     if (verbose) PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",rank,__FUNCT__);
766     ierr = PetscFree(clid_flid);CHKERRQ(ierr);
767     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
768     *a_P_out = NULL;  /* out */
769     PetscFunctionReturn(0);
770   }
771 
772   {
773     PetscReal *coords;
774     PetscInt  data_stride;
775     PetscInt  *crsGID = NULL;
776     Mat       Gmat2;
777 
778     if (dim != data_cols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim %D != data_cols %D",dim,data_cols);
779     /* grow ghost data for better coarse grid cover of fine grid */
780 #if defined PETSC_GAMG_USE_LOG
781     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
782 #endif
783     /* messy method, squares graph and gets some data */
784     ierr = getGIDsOnSquareGraph(nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID);CHKERRQ(ierr);
785     /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */
786 #if defined PETSC_GAMG_USE_LOG
787     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr);
788 #endif
789     /* create global vector of coorindates in 'coords' */
790     if (size > 1) {
791       ierr = PCGAMGGetDataWithGhosts(Gmat2, dim, pc_gamg->data, &data_stride, &coords);CHKERRQ(ierr);
792     } else {
793       coords      = (PetscReal*)pc_gamg->data;
794       data_stride = pc_gamg->data_sz/pc_gamg->data_cell_cols;
795     }
796     ierr = MatDestroy(&Gmat2);CHKERRQ(ierr);
797 
798     /* triangulate */
799     if (dim == 2) {
800       PetscReal metric,tm;
801 #if defined PETSC_GAMG_USE_LOG
802       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
803 #endif
804       ierr = triangulateAndFormProl(selected_2, data_stride, coords,nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric);CHKERRQ(ierr);
805 #if defined PETSC_GAMG_USE_LOG
806       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr);
807 #endif
808       ierr = PetscFree(crsGID);CHKERRQ(ierr);
809 
810       /* clean up and create coordinates for coarse grid (output) */
811       if (size > 1) ierr = PetscFree(coords);CHKERRQ(ierr);
812 
813       ierr = MPI_Allreduce(&metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm);CHKERRQ(ierr);
814       if (tm > 1.) { /* needs to be globalized - should not happen */
815         if (verbose) PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",rank,__FUNCT__,tm);
816         ierr = MatDestroy(&Prol);CHKERRQ(ierr);
817         Prol = NULL;
818       } else if (metric > .0) {
819         if (verbose) PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",rank,__FUNCT__,metric);
820       }
821     } else SETERRQ(wcomm,PETSC_ERR_PLIB,"3D not implemented for 'geo' AMG");
822     { /* create next coords - output */
823       PetscReal *crs_crds;
824       ierr = PetscMalloc(dim*nLocalSelected*sizeof(PetscReal), &crs_crds);CHKERRQ(ierr);
825       for (kk=0; kk<nLocalSelected; kk++) { /* grab local select nodes to promote - output */
826         PetscInt lid = clid_flid[kk];
827         for (jj=0; jj<dim; jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid];
828       }
829 
830       ierr             = PetscFree(pc_gamg->data);CHKERRQ(ierr);
831       pc_gamg->data    = crs_crds; /* out */
832       pc_gamg->data_sz = dim*nLocalSelected;
833     }
834     ierr = ISDestroy(&selected_2);CHKERRQ(ierr);
835   }
836 
837   *a_P_out = Prol;  /* out */
838   ierr     = PetscFree(clid_flid);CHKERRQ(ierr);
839 #if defined PETSC_USE_LOG
840   ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr);
841 #endif
842   PetscFunctionReturn(0);
843 }
844 
845 /* -------------------------------------------------------------------------- */
846 /*
847  PCCreateGAMG_GEO
848 
849   Input Parameter:
850    . pc -
851 */
852 #undef __FUNCT__
853 #define __FUNCT__ "PCCreateGAMG_GEO"
854 PetscErrorCode  PCCreateGAMG_GEO(PC pc)
855 {
856   PetscErrorCode ierr;
857   PC_MG          *mg      = (PC_MG*)pc->data;
858   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
859 
860   PetscFunctionBegin;
861   pc->ops->setfromoptions = PCSetFromOptions_GEO;
862   /* pc->ops->destroy        = PCDestroy_GEO; */
863   /* reset does not do anything; setup not virtual */
864 
865   /* set internal function pointers */
866   pc_gamg->graph       = PCGAMGgraph_GEO;
867   pc_gamg->coarsen     = PCGAMGcoarsen_GEO;
868   pc_gamg->prolongator = PCGAMGProlongator_GEO;
869   pc_gamg->optprol     = 0;
870   pc_gamg->formkktprol = 0;
871 
872   pc_gamg->createdefaultdata = PCSetData_GEO;
873 
874   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,
875                                            "PCSetCoordinates_C",
876                                            "PCSetCoordinates_GEO",
877                                            PCSetCoordinates_GEO);CHKERRQ(ierr);
878   PetscFunctionReturn(0);
879 }
880