Actual source code: bddcgraph.c

  1: #include <petsc/private/petscimpl.h>
  2: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
  3: #include <../src/ksp/pc/impls/bddc/bddcstructs.h>

  5: PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS* dirdofs)
  6: {
  7:   if (graph->dirdofsB) {
  8:     PetscObjectReference((PetscObject)graph->dirdofsB);
  9:   } else if (graph->has_dirichlet) {
 10:     PetscInt i,size;
 11:     PetscInt *dirdofs_idxs;

 13:     size = 0;
 14:     for (i=0;i<graph->nvtxs;i++) {
 15:       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
 16:     }

 18:     PetscMalloc1(size,&dirdofs_idxs);
 19:     size = 0;
 20:     for (i=0;i<graph->nvtxs;i++) {
 21:       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
 22:     }
 23:     ISCreateGeneral(PETSC_COMM_SELF,size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofsB);
 24:     PetscObjectReference((PetscObject)graph->dirdofsB);
 25:   }
 26:   *dirdofs = graph->dirdofsB;
 27:   return 0;
 28: }

 30: PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
 31: {
 32:   if (graph->dirdofs) {
 33:     PetscObjectReference((PetscObject)graph->dirdofs);
 34:   } else if (graph->has_dirichlet) {
 35:     PetscInt i,size;
 36:     PetscInt *dirdofs_idxs;

 38:     size = 0;
 39:     for (i=0;i<graph->nvtxs;i++) {
 40:       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
 41:     }

 43:     PetscMalloc1(size,&dirdofs_idxs);
 44:     size = 0;
 45:     for (i=0;i<graph->nvtxs;i++) {
 46:       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
 47:     }
 48:     ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs);
 49:     PetscObjectReference((PetscObject)graph->dirdofs);
 50:   }
 51:   *dirdofs = graph->dirdofs;
 52:   return 0;
 53: }

 55: PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
 56: {
 57:   PetscInt       i,j,tabs;
 58:   PetscInt*      queue_in_global_numbering;

 60:   PetscViewerASCIIPushSynchronized(viewer);
 61:   PetscViewerASCIIGetTab(viewer,&tabs);
 62:   PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");
 63:   PetscViewerFlush(viewer);
 64:   PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);
 65:   PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);
 66:   PetscViewerASCIISynchronizedPrintf(viewer,"Number of local subdomains %d\n",graph->n_local_subs ? graph->n_local_subs : 1);
 67:   PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);
 68:   if (graph->maxcount != PETSC_MAX_INT) {
 69:     PetscViewerASCIISynchronizedPrintf(viewer,"Max count %d\n",graph->maxcount);
 70:   }
 71:   PetscViewerASCIISynchronizedPrintf(viewer,"Topological two dim? %d (set %d)\n",graph->twodim,graph->twodimset);
 72:   if (verbosity_level > 2) {
 73:     for (i=0;i<graph->nvtxs;i++) {
 74:       PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);
 75:       PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);
 76:       PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);
 77:       PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);
 78:       PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
 79:       if (graph->count[i]) {
 80:         PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");
 81:         for (j=0;j<graph->count[i];j++) {
 82:           PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);
 83:         }
 84:         PetscViewerASCIISynchronizedPrintf(viewer,"\n");
 85:       }
 86:       PetscViewerASCIISetTab(viewer,tabs);
 87:       PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
 88:       if (graph->mirrors) {
 89:         PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);
 90:         if (graph->mirrors[i]) {
 91:           PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
 92:           PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");
 93:           for (j=0;j<graph->mirrors[i];j++) {
 94:             PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);
 95:           }
 96:           PetscViewerASCIISynchronizedPrintf(viewer,"\n");
 97:           PetscViewerASCIISetTab(viewer,tabs);
 98:           PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
 99:         }
100:       }
101:       if (verbosity_level > 3) {
102:         if (graph->xadj) {
103:           PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");
104:           PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
105:           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
106:             PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);
107:           }
108:           PetscViewerASCIISynchronizedPrintf(viewer,"\n");
109:           PetscViewerASCIISetTab(viewer,tabs);
110:           PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
111:         } else {
112:           PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n");
113:         }
114:       }
115:       if (graph->n_local_subs) {
116:         PetscViewerASCIISynchronizedPrintf(viewer,"   local sub id: %d\n",graph->local_subs[i]);
117:       }
118:       PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);
119:       if (graph->subset[i] && graph->subset_ncc) {
120:         PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);
121:       }
122:     }
123:   }
124:   PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);
125:   PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);
126:   ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);
127:   for (i=0;i<graph->ncc;i++) {
128:     PetscInt node_num=graph->queue[graph->cptr[i]];
129:     PetscBool printcc = PETSC_FALSE;
130:     PetscViewerASCIISynchronizedPrintf(viewer,"  cc %d (size %d, fid %d, neighs:",i,graph->cptr[i+1]-graph->cptr[i],graph->which_dof[node_num]);
131:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
132:     for (j=0;j<graph->count[node_num];j++) {
133:       PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);
134:     }
135:     if (verbosity_level > 1) {
136:       PetscViewerASCIISynchronizedPrintf(viewer,"):");
137:       if (verbosity_level > 2 || graph->twodim || graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
138:         printcc = PETSC_TRUE;
139:       }
140:       if (printcc) {
141:         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
142:           PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);
143:         }
144:       }
145:     } else {
146:       PetscViewerASCIISynchronizedPrintf(viewer,")");
147:     }
148:     PetscViewerASCIISynchronizedPrintf(viewer,"\n");
149:     PetscViewerASCIISetTab(viewer,tabs);
150:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
151:   }
152:   PetscFree(queue_in_global_numbering);
153:   PetscViewerFlush(viewer);
154:   return 0;
155: }

157: PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
158: {
159:   PetscInt       i;

161:   if (n_faces) {
162:     if (FacesIS) {
163:       for (i=0;i<*n_faces;i++) {
164:         ISDestroy(&((*FacesIS)[i]));
165:       }
166:       PetscFree(*FacesIS);
167:     }
168:     *n_faces = 0;
169:   }
170:   if (n_edges) {
171:     if (EdgesIS) {
172:       for (i=0;i<*n_edges;i++) {
173:         ISDestroy(&((*EdgesIS)[i]));
174:       }
175:       PetscFree(*EdgesIS);
176:     }
177:     *n_edges = 0;
178:   }
179:   if (VerticesIS) {
180:     ISDestroy(VerticesIS);
181:   }
182:   return 0;
183: }

185: PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
186: {
187:   IS             *ISForFaces,*ISForEdges,ISForVertices;
188:   PetscInt       i,nfc,nec,nvc,*idx,*mark;

190:   PetscCalloc1(graph->ncc,&mark);
191:   /* loop on ccs to evalute number of faces, edges and vertices */
192:   nfc = 0;
193:   nec = 0;
194:   nvc = 0;
195:   for (i=0;i<graph->ncc;i++) {
196:     PetscInt repdof = graph->queue[graph->cptr[i]];
197:     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size && graph->count[repdof] < graph->maxcount) {
198:       if (!graph->twodim && graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
199:         nfc++;
200:         mark[i] = 2;
201:       } else {
202:         nec++;
203:         mark[i] = 1;
204:       }
205:     } else {
206:       nvc += graph->cptr[i+1]-graph->cptr[i];
207:     }
208:   }

210:   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
211:   if (FacesIS) {
212:     PetscMalloc1(nfc,&ISForFaces);
213:   }
214:   if (EdgesIS) {
215:     PetscMalloc1(nec,&ISForEdges);
216:   }
217:   if (VerticesIS) {
218:     PetscMalloc1(nvc,&idx);
219:   }

221:   /* loop on ccs to compute index sets for faces and edges */
222:   if (!graph->queue_sorted) {
223:     PetscInt *queue_global;

225:     PetscMalloc1(graph->cptr[graph->ncc],&queue_global);
226:     ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);
227:     for (i=0;i<graph->ncc;i++) {
228:       PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);
229:     }
230:     PetscFree(queue_global);
231:     graph->queue_sorted = PETSC_TRUE;
232:   }
233:   nfc = 0;
234:   nec = 0;
235:   for (i=0;i<graph->ncc;i++) {
236:     if (mark[i] == 2) {
237:       if (FacesIS) {
238:         ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForFaces[nfc]);
239:       }
240:       nfc++;
241:     } else if (mark[i] == 1) {
242:       if (EdgesIS) {
243:         ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForEdges[nec]);
244:       }
245:       nec++;
246:     }
247:   }

249:   /* index set for vertices */
250:   if (VerticesIS) {
251:     nvc = 0;
252:     for (i=0;i<graph->ncc;i++) {
253:       if (!mark[i]) {
254:         PetscInt j;

256:         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
257:           idx[nvc]=graph->queue[j];
258:           nvc++;
259:         }
260:       }
261:     }
262:     /* sort vertex set (by local ordering) */
263:     PetscSortInt(nvc,idx);
264:     ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);
265:   }
266:   PetscFree(mark);

268:   /* get back info */
269:   if (n_faces)       *n_faces = nfc;
270:   if (FacesIS)       *FacesIS = ISForFaces;
271:   if (n_edges)       *n_edges = nec;
272:   if (EdgesIS)       *EdgesIS = ISForEdges;
273:   if (VerticesIS) *VerticesIS = ISForVertices;
274:   return 0;
275: }

277: PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
278: {
279:   PetscBool      adapt_interface_reduced;
280:   MPI_Comm       interface_comm;
281:   PetscMPIInt    size;
282:   PetscInt       i;
283:   PetscBT        cornerp;

285:   /* compute connected components locally */
286:   PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);
287:   PCBDDCGraphComputeConnectedComponentsLocal(graph);

289:   cornerp = NULL;
290:   if (graph->active_coords) { /* face based corner selection */
291:     PetscBT   excluded;
292:     PetscReal *wdist;
293:     PetscInt  n_neigh,*neigh,*n_shared,**shared;
294:     PetscInt  maxc, ns;

296:     PetscBTCreate(graph->nvtxs,&cornerp);
297:     ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);
298:     for (ns = 1, maxc = 0; ns < n_neigh; ns++) maxc = PetscMax(maxc,n_shared[ns]);
299:     PetscMalloc1(maxc*graph->cdim,&wdist);
300:     PetscBTCreate(maxc,&excluded);

302:     for (ns = 1; ns < n_neigh; ns++) { /* first proc is self */
303:       PetscReal *anchor,mdist;
304:       PetscInt  fst,j,k,d,cdim = graph->cdim,n = n_shared[ns];
305:       PetscInt  point1,point2,point3,point4;

307:       /* import coordinates on shared interface */
308:       PetscBTMemzero(n,excluded);
309:       for (j=0,fst=-1,k=0;j<n;j++) {
310:         PetscBool skip = PETSC_FALSE;
311:         for (d=0;d<cdim;d++) {
312:           PetscReal c = graph->coords[shared[ns][j]*cdim+d];
313:           skip = (PetscBool)(skip || c == PETSC_MAX_REAL);
314:           wdist[k++] = c;
315:         }
316:         if (skip) {
317:           PetscBTSet(excluded,j);
318:         } else if (fst == -1) fst = j;
319:       }
320:       if (fst == -1) continue;

322:       /* the dofs are sorted by global numbering, so each rank starts from the same id
323:          and it will detect the same corners from the given set */

325:       /* find the farthest point from the starting one */
326:       anchor = wdist + fst*cdim;
327:       mdist  = -1.0;
328:       point1 = fst;
329:       for (j=fst;j<n;j++) {
330:         PetscReal dist = 0.0;

332:         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
333:         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
334:         if (dist > mdist) { mdist = dist; point1 = j; }
335:       }

337:       /* find the farthest point from point1 */
338:       anchor = wdist + point1*cdim;
339:       mdist  = -1.0;
340:       point2 = point1;
341:       for (j=fst;j<n;j++) {
342:         PetscReal dist = 0.0;

344:         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
345:         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
346:         if (dist > mdist) { mdist = dist; point2 = j; }
347:       }

349:       /* find the third point maximizing the triangle area */
350:       point3 = point2;
351:       if (cdim > 2) {
352:         PetscReal a = 0.0;

354:         for (d=0;d<cdim;d++) a += (wdist[point1*cdim+d]-wdist[point2*cdim+d])*(wdist[point1*cdim+d]-wdist[point2*cdim+d]);
355:         a = PetscSqrtReal(a);
356:         mdist = -1.0;
357:         for (j=fst;j<n;j++) {
358:           PetscReal area,b = 0.0, c = 0.0,s;

360:           if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
361:           for (d=0;d<cdim;d++) {
362:             b += (wdist[point1*cdim+d]-wdist[j*cdim+d])*(wdist[point1*cdim+d]-wdist[j*cdim+d]);
363:             c += (wdist[point2*cdim+d]-wdist[j*cdim+d])*(wdist[point2*cdim+d]-wdist[j*cdim+d]);
364:           }
365:           b = PetscSqrtReal(b);
366:           c = PetscSqrtReal(c);
367:           s = 0.5*(a+b+c);

369:           /* Heron's formula, area squared */
370:           area = s*(s-a)*(s-b)*(s-c);
371:           if (area > mdist) { mdist = area; point3 = j; }
372:         }
373:       }

375:       /* find the farthest point from point3 different from point1 and point2 */
376:       anchor = wdist + point3*cdim;
377:       mdist  = -1.0;
378:       point4 = point3;
379:       for (j=fst;j<n;j++) {
380:         PetscReal dist = 0.0;

382:         if (PetscUnlikely(PetscBTLookup(excluded,j)) || j == point1 || j == point2 || j == point3) continue;
383:         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
384:         if (dist > mdist) { mdist = dist; point4 = j; }
385:       }

387:       PetscBTSet(cornerp,shared[ns][point1]);
388:       PetscBTSet(cornerp,shared[ns][point2]);
389:       PetscBTSet(cornerp,shared[ns][point3]);
390:       PetscBTSet(cornerp,shared[ns][point4]);

392:       /* all dofs having the same coordinates will be primal */
393:       for (j=fst;j<n;j++) {
394:         PetscBool same[] = {PETSC_TRUE,PETSC_TRUE,PETSC_TRUE,PETSC_TRUE};

396:         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
397:         for (d=0;d<cdim;d++) {
398:           same[0] = (PetscBool)(same[0] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point1*cdim+d]) < PETSC_SMALL));
399:           same[1] = (PetscBool)(same[1] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point2*cdim+d]) < PETSC_SMALL));
400:           same[2] = (PetscBool)(same[2] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point3*cdim+d]) < PETSC_SMALL));
401:           same[3] = (PetscBool)(same[3] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point4*cdim+d]) < PETSC_SMALL));
402:         }
403:         if (same[0] || same[1] || same[2] || same[3]) {
404:           PetscBTSet(cornerp,shared[ns][j]);
405:         }
406:       }
407:     }
408:     PetscBTDestroy(&excluded);
409:     PetscFree(wdist);
410:     ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);
411:   }

413:   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
414:   MPI_Comm_size(interface_comm,&size);
415:   adapt_interface_reduced = PETSC_FALSE;
416:   if (size > 1) {
417:     PetscInt i;
418:     PetscBool adapt_interface = cornerp ? PETSC_TRUE : PETSC_FALSE;
419:     for (i=0;i<graph->n_subsets && !adapt_interface;i++) {
420:       /* We are not sure that on a given subset of the local interface,
421:          with two connected components, the latters be the same among sharing subdomains */
422:       if (graph->subset_ncc[i] > 1) adapt_interface = PETSC_TRUE;
423:     }
424:     MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);
425:   }

427:   if (graph->n_subsets && adapt_interface_reduced) {
428:     PetscBT     subset_cc_adapt;
429:     MPI_Request *send_requests,*recv_requests;
430:     PetscInt    *send_buffer,*recv_buffer;
431:     PetscInt    sum_requests,start_of_recv,start_of_send;
432:     PetscInt    *cum_recv_counts;
433:     PetscInt    *labels;
434:     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
435:     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
436:     PetscBool   *subset_has_corn,*recv_buffer_bool,*send_buffer_bool;

438:     PetscCalloc1(graph->n_subsets,&subset_has_corn);
439:     if (cornerp) {
440:       for (i=0;i<graph->n_subsets;i++) {
441:         for (j=0;j<graph->subset_size[i];j++) {
442:           if (PetscBTLookup(cornerp,graph->subset_idxs[i][j])) {
443:             subset_has_corn[i] = PETSC_TRUE;
444:             break;
445:           }
446:         }
447:       }
448:     }
449:     PetscMalloc1(graph->nvtxs,&labels);
450:     PetscArrayzero(labels,graph->nvtxs);
451:     for (i=0,k=0;i<graph->ncc;i++) {
452:       PetscInt s = 1;
453:       for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
454:         if (cornerp && PetscBTLookup(cornerp,graph->queue[j])) {
455:           labels[graph->queue[j]] = k+s;
456:           s += 1;
457:         } else {
458:           labels[graph->queue[j]] = k;
459:         }
460:       }
461:       k += s;
462:     }

464:     /* allocate some space */
465:     PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);
466:     PetscArrayzero(cum_recv_counts,graph->n_subsets+1);

468:     /* first count how many neighbours per connected component I will receive from */
469:     cum_recv_counts[0] = 0;
470:     for (i=0;i<graph->n_subsets;i++) cum_recv_counts[i+1] = cum_recv_counts[i]+graph->count[graph->subset_idxs[i][0]];
471:     PetscMalloc1(graph->n_subsets,&send_buffer_bool);
472:     PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer_bool);
473:     PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);
474:     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
475:       send_requests[i] = MPI_REQUEST_NULL;
476:       recv_requests[i] = MPI_REQUEST_NULL;
477:     }

479:     /* exchange with my neighbours the number of my connected components on the subset of interface */
480:     sum_requests = 0;
481:     for (i=0;i<graph->n_subsets;i++) {
482:       send_buffer_bool[i] = (PetscBool)(graph->subset_ncc[i] > 1 || subset_has_corn[i]);
483:     }
484:     for (i=0;i<graph->n_subsets;i++) {
485:       PetscMPIInt neigh,tag;
486:       PetscInt    count,*neighs;

488:       count  = graph->count[graph->subset_idxs[i][0]];
489:       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
490:       PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);
491:       for (k=0;k<count;k++) {

493:         PetscMPIIntCast(neighs[k],&neigh);
494:         MPI_Isend(send_buffer_bool + i,           1,MPIU_BOOL,neigh,tag,interface_comm,&send_requests[sum_requests]);
495:         MPI_Irecv(recv_buffer_bool + sum_requests,1,MPIU_BOOL,neigh,tag,interface_comm,&recv_requests[sum_requests]);
496:         sum_requests++;
497:       }
498:     }
499:     MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);
500:     MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);

502:     /* determine the subsets I have to adapt (those having more than 1 cc) */
503:     PetscBTCreate(graph->n_subsets,&subset_cc_adapt);
504:     PetscBTMemzero(graph->n_subsets,subset_cc_adapt);
505:     for (i=0;i<graph->n_subsets;i++) {
506:       if (graph->subset_ncc[i] > 1 || subset_has_corn[i]) {
507:         PetscBTSet(subset_cc_adapt,i);
508:         continue;
509:       }
510:       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++) {
511:          if (recv_buffer_bool[j]) {
512:           PetscBTSet(subset_cc_adapt,i);
513:           break;
514:         }
515:       }
516:     }
517:     PetscFree(send_buffer_bool);
518:     PetscFree(recv_buffer_bool);
519:     PetscFree(subset_has_corn);

521:     /* determine send/recv buffers sizes */
522:     j = 0;
523:     mss = 0;
524:     for (i=0;i<graph->n_subsets;i++) {
525:       if (PetscBTLookup(subset_cc_adapt,i)) {
526:         j  += graph->subset_size[i];
527:         mss = PetscMax(graph->subset_size[i],mss);
528:       }
529:     }
530:     k = 0;
531:     mns = 0;
532:     for (i=0;i<graph->n_subsets;i++) {
533:       if (PetscBTLookup(subset_cc_adapt,i)) {
534:         k  += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
535:         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
536:       }
537:     }
538:     PetscMalloc2(j,&send_buffer,k,&recv_buffer);

540:     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
541:     j = 0;
542:     for (i=0;i<graph->n_subsets;i++)
543:       if (PetscBTLookup(subset_cc_adapt,i))
544:         for (k=0;k<graph->subset_size[i];k++)
545:           send_buffer[j++] = labels[graph->subset_idxs[i][k]];

547:     /* now exchange the data */
548:     start_of_recv = 0;
549:     start_of_send = 0;
550:     sum_requests  = 0;
551:     for (i=0;i<graph->n_subsets;i++) {
552:       if (PetscBTLookup(subset_cc_adapt,i)) {
553:         PetscMPIInt neigh,tag;
554:         PetscInt    size_of_send = graph->subset_size[i];

556:         j    = graph->subset_idxs[i][0];
557:         PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);
558:         for (k=0;k<graph->count[j];k++) {
559:           PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);
560:           MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);
561:           MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);
562:           start_of_recv += size_of_send;
563:           sum_requests++;
564:         }
565:         start_of_send += size_of_send;
566:       }
567:     }
568:     MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);

570:     /* refine connected components */
571:     start_of_recv = 0;
572:     /* allocate some temporary space */
573:     if (mss) {
574:       PetscMalloc1(mss,&refine_buffer);
575:       PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);
576:     }
577:     ncc = 0;
578:     cum_queue = 0;
579:     graph->cptr[0] = 0;
580:     for (i=0;i<graph->n_subsets;i++) {
581:       if (PetscBTLookup(subset_cc_adapt,i)) {
582:         PetscInt subset_counter = 0;
583:         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
584:         PetscInt buffer_size = graph->subset_size[i];

586:         /* compute pointers */
587:         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
588:         /* analyze contributions from subdomains that share the i-th subset
589:            The structure of refine_buffer is suitable to find intersections of ccs among sharingprocs.
590:            supposing the current subset is shared by 3 processes and has dimension 5 with global dofs 0,1,2,3,4 (local 0,4,3,1,2)
591:            sharing procs connected components:
592:              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
593:              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
594:              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
595:            refine_buffer will be filled as:
596:              [ 4, 3, 1;
597:                4, 2, 1;
598:                7, 2, 6;
599:                4, 3, 5;
600:                7, 2, 6; ];
601:            The connected components in local ordering are [0], [1], [2 3], [4] */
602:         /* fill temp_buffer */
603:         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
604:         for (j=0;j<sharingprocs-1;j++) {
605:           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
606:           start_of_recv += buffer_size;
607:         }
608:         PetscArrayzero(private_labels,buffer_size);
609:         for (j=0;j<buffer_size;j++) {
610:           if (!private_labels[j]) { /* found a new cc  */
611:             PetscBool same_set;

613:             graph->cptr[ncc] = cum_queue;
614:             ncc++;
615:             subset_counter++;
616:             private_labels[j] = subset_counter;
617:             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
618:             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
619:               same_set = PETSC_TRUE;
620:               for (s=0;s<sharingprocs;s++) {
621:                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
622:                   same_set = PETSC_FALSE;
623:                   break;
624:                 }
625:               }
626:               if (same_set) {
627:                 private_labels[k] = subset_counter;
628:                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
629:               }
630:             }
631:           }
632:         }
633:         graph->cptr[ncc]     = cum_queue;
634:         graph->subset_ncc[i] = subset_counter;
635:         graph->queue_sorted  = PETSC_FALSE;
636:       } else { /* this subset does not need to be adapted */
637:         PetscArraycpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]);
638:         ncc++;
639:         cum_queue += graph->subset_size[i];
640:         graph->cptr[ncc] = cum_queue;
641:       }
642:     }
643:     graph->cptr[ncc] = cum_queue;
644:     graph->ncc       = ncc;
645:     if (mss) {
646:       PetscFree2(refine_buffer[0],private_labels);
647:       PetscFree(refine_buffer);
648:     }
649:     PetscFree(labels);
650:     MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);
651:     PetscFree2(send_requests,recv_requests);
652:     PetscFree2(send_buffer,recv_buffer);
653:     PetscFree(cum_recv_counts);
654:     PetscBTDestroy(&subset_cc_adapt);
655:   }
656:   PetscBTDestroy(&cornerp);

658:   /* Determine if we are in 2D or 3D */
659:   if (!graph->twodimset) {
660:     PetscBool twodim = PETSC_TRUE;
661:     for (i=0;i<graph->ncc;i++) {
662:       PetscInt repdof = graph->queue[graph->cptr[i]];
663:       PetscInt ccsize = graph->cptr[i+1]-graph->cptr[i];
664:       if (graph->count[repdof] > 1 && ccsize > graph->custom_minimal_size) {
665:         twodim = PETSC_FALSE;
666:         break;
667:       }
668:     }
669:     MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));
670:     graph->twodimset = PETSC_TRUE;
671:   }
672:   return 0;
673: }

675: static inline PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
676: {
677:   PetscInt       i,j,n;
678:   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
679:   PetscBT        touched = graph->touched;
680:   PetscBool      havecsr = (PetscBool)(!!xadj);
681:   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);

683:   n = 0;
684:   if (havecsr && !havesubs) {
685:     for (i=-n_prev;i<0;i++) {
686:       PetscInt start_dof = queue_tip[i];
687:       /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs */
688:       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
689:         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
690:           PetscInt dof = graph->subset_idxs[pid-1][j];
691:           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
692:             PetscBTSet(touched,dof);
693:             queue_tip[n] = dof;
694:             n++;
695:           }
696:         }
697:       } else {
698:         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
699:           PetscInt dof = adjncy[j];
700:           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
701:             PetscBTSet(touched,dof);
702:             queue_tip[n] = dof;
703:             n++;
704:           }
705:         }
706:       }
707:     }
708:   } else if (havecsr && havesubs) {
709:     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
710:     for (i=-n_prev;i<0;i++) {
711:       PetscInt start_dof = queue_tip[i];
712:       /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs belonging to the local sub */
713:       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
714:         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
715:           PetscInt dof = graph->subset_idxs[pid-1][j];
716:           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
717:             PetscBTSet(touched,dof);
718:             queue_tip[n] = dof;
719:             n++;
720:           }
721:         }
722:       } else {
723:         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
724:           PetscInt dof = adjncy[j];
725:           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
726:             PetscBTSet(touched,dof);
727:             queue_tip[n] = dof;
728:             n++;
729:           }
730:         }
731:       }
732:     }
733:   } else if (havesubs) { /* sub info only */
734:     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
735:     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
736:       PetscInt dof = graph->subset_idxs[pid-1][j];
737:       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
738:         PetscBTSet(touched,dof);
739:         queue_tip[n] = dof;
740:         n++;
741:       }
742:     }
743:   } else {
744:     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
745:       PetscInt dof = graph->subset_idxs[pid-1][j];
746:       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
747:         PetscBTSet(touched,dof);
748:         queue_tip[n] = dof;
749:         n++;
750:       }
751:     }
752:   }
753:   *n_added = n;
754:   return 0;
755: }

757: PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
758: {
759:   PetscInt       ncc,cum_queue,n;
760:   PetscMPIInt    commsize;

763:   /* quiet return if there isn't any local info */
764:   if (!graph->xadj && !graph->n_local_subs) {
765:     return 0;
766:   }

768:   /* reset any previous search of connected components */
769:   PetscBTMemzero(graph->nvtxs,graph->touched);
770:   MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);
771:   if (commsize > graph->commsizelimit) {
772:     PetscInt i;
773:     for (i=0;i<graph->nvtxs;i++) {
774:       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
775:         PetscBTSet(graph->touched,i);
776:       }
777:     }
778:   }

780:   /* begin search for connected components */
781:   cum_queue = 0;
782:   ncc = 0;
783:   for (n=0;n<graph->n_subsets;n++) {
784:     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
785:     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
786:     while (found != graph->subset_size[n]) {
787:       PetscInt added = 0;
788:       if (!prev) { /* search for new starting dof */
789:         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
790:         PetscBTSet(graph->touched,graph->subset_idxs[n][first]);
791:         graph->queue[cum_queue] = graph->subset_idxs[n][first];
792:         graph->cptr[ncc] = cum_queue;
793:         prev = 1;
794:         cum_queue++;
795:         found++;
796:         ncc_pid++;
797:         ncc++;
798:       }
799:       PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);
800:       if (!added) {
801:         graph->subset_ncc[n] = ncc_pid;
802:         graph->cptr[ncc] = cum_queue;
803:       }
804:       prev = added;
805:       found += added;
806:       cum_queue += added;
807:       if (added && found == graph->subset_size[n]) {
808:         graph->subset_ncc[n] = ncc_pid;
809:         graph->cptr[ncc] = cum_queue;
810:       }
811:     }
812:   }
813:   graph->ncc = ncc;
814:   graph->queue_sorted = PETSC_FALSE;
815:   return 0;
816: }

818: PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
819: {
820:   IS             subset,subset_n;
821:   MPI_Comm       comm;
822:   const PetscInt *is_indices;
823:   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
824:   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
825:   PetscMPIInt    commsize;
826:   PetscBool      same_set,mirrors_found;

829:   if (neumann_is) {
832:   }
833:   graph->has_dirichlet = PETSC_FALSE;
834:   if (dirichlet_is) {
837:     graph->has_dirichlet = PETSC_TRUE;
838:   }
840:   for (i=0;i<n_ISForDofs;i++) {
843:   }
844:   if (custom_primal_vertices) {
847:   }
848:   PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);
849:   MPI_Comm_size(comm,&commsize);

851:   /* custom_minimal_size */
852:   graph->custom_minimal_size = custom_minimal_size;
853:   /* get info l2gmap and allocate work vectors  */
854:   ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);
855:   /* check if we have any local periodic nodes (periodic BCs) */
856:   mirrors_found = PETSC_FALSE;
857:   if (graph->nvtxs && n_neigh) {
858:     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
859:     for (i=0; i<n_shared[0]; i++) {
860:       if (graph->count[shared[0][i]] > 1) {
861:         mirrors_found = PETSC_TRUE;
862:         break;
863:       }
864:     }
865:   }
866:   /* compute local mirrors (if any) */
867:   if (mirrors_found) {
868:     IS       to,from;
869:     PetscInt *local_indices,*global_indices;

871:     ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);
872:     ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);
873:     /* get arrays of local and global indices */
874:     PetscMalloc1(graph->nvtxs,&local_indices);
875:     ISGetIndices(to,(const PetscInt**)&is_indices);
876:     PetscArraycpy(local_indices,is_indices,graph->nvtxs);
877:     ISRestoreIndices(to,(const PetscInt**)&is_indices);
878:     PetscMalloc1(graph->nvtxs,&global_indices);
879:     ISGetIndices(from,(const PetscInt**)&is_indices);
880:     PetscArraycpy(global_indices,is_indices,graph->nvtxs);
881:     ISRestoreIndices(from,(const PetscInt**)&is_indices);
882:     /* allocate space for mirrors */
883:     PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);
884:     PetscArrayzero(graph->mirrors,graph->nvtxs);
885:     graph->mirrors_set[0] = NULL;

887:     k=0;
888:     for (i=0;i<n_shared[0];i++) {
889:       j=shared[0][i];
890:       if (graph->count[j] > 1) {
891:         graph->mirrors[j]++;
892:         k++;
893:       }
894:     }
895:     /* allocate space for set of mirrors */
896:     PetscMalloc1(k,&graph->mirrors_set[0]);
897:     for (i=1;i<graph->nvtxs;i++)
898:       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];

900:     /* fill arrays */
901:     PetscArrayzero(graph->mirrors,graph->nvtxs);
902:     for (j=0;j<n_shared[0];j++) {
903:       i=shared[0][j];
904:       if (graph->count[i] > 1)
905:         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
906:     }
907:     PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);
908:     for (i=0;i<graph->nvtxs;i++) {
909:       if (graph->mirrors[i] > 0) {
910:         PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);
911:         j = global_indices[k];
912:         while (k > 0 && global_indices[k-1] == j) k--;
913:         for (j=0;j<graph->mirrors[i];j++) {
914:           graph->mirrors_set[i][j]=local_indices[k+j];
915:         }
916:         PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);
917:       }
918:     }
919:     PetscFree(local_indices);
920:     PetscFree(global_indices);
921:     ISDestroy(&to);
922:     ISDestroy(&from);
923:   }
924:   PetscArrayzero(graph->count,graph->nvtxs);

926:   /* Count total number of neigh per node */
927:   k = 0;
928:   for (i=1;i<n_neigh;i++) {
929:     k += n_shared[i];
930:     for (j=0;j<n_shared[i];j++) {
931:       graph->count[shared[i][j]] += 1;
932:     }
933:   }
934:   /* Allocate space for storing the set of neighbours for each node */
935:   if (graph->nvtxs) {
936:     PetscMalloc1(k,&graph->neighbours_set[0]);
937:   }
938:   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
939:     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
940:   }
941:   /* Get information for sharing subdomains */
942:   PetscArrayzero(graph->count,graph->nvtxs);
943:   for (i=1;i<n_neigh;i++) { /* dont count myself */
944:     s = n_shared[i];
945:     for (j=0;j<s;j++) {
946:       k = shared[i][j];
947:       graph->neighbours_set[k][graph->count[k]] = neigh[i];
948:       graph->count[k] += 1;
949:     }
950:   }
951:   /* sort set of sharing subdomains */
952:   for (i=0;i<graph->nvtxs;i++) {
953:     PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);
954:   }
955:   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
956:   ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);

958:   /*
959:      Get info for dofs splitting
960:      User can specify just a subset; an additional field is considered as a complementary field
961:   */
962:   for (i=0,k=0;i<n_ISForDofs;i++) {
963:     PetscInt bs;

965:     ISGetBlockSize(ISForDofs[i],&bs);
966:     k   += bs;
967:   }
968:   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = k; /* by default a dof belongs to the complement set */
969:   for (i=0,k=0;i<n_ISForDofs;i++) {
970:     PetscInt bs;

972:     ISGetLocalSize(ISForDofs[i],&is_size);
973:     ISGetBlockSize(ISForDofs[i],&bs);
974:     ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);
975:     for (j=0;j<is_size/bs;j++) {
976:       PetscInt b;

978:       for (b=0;b<bs;b++) {
979:         PetscInt jj = bs*j + b;

981:         if (is_indices[jj] > -1 && is_indices[jj] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
982:           graph->which_dof[is_indices[jj]] = k+b;
983:         }
984:       }
985:     }
986:     ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);
987:     k   += bs;
988:   }

990:   /* Take into account Neumann nodes */
991:   if (neumann_is) {
992:     ISGetLocalSize(neumann_is,&is_size);
993:     ISGetIndices(neumann_is,(const PetscInt**)&is_indices);
994:     for (i=0;i<is_size;i++) {
995:       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
996:         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_NEUMANN_MARK;
997:       }
998:     }
999:     ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);
1000:   }
1001:   /* Take into account Dirichlet nodes (they overwrite any neumann boundary mark previously set) */
1002:   if (dirichlet_is) {
1003:     ISGetLocalSize(dirichlet_is,&is_size);
1004:     ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);
1005:     for (i=0;i<is_size;i++) {
1006:       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
1007:         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
1008:           PetscBTSet(graph->touched,is_indices[i]);
1009:           graph->subset[is_indices[i]] = 0;
1010:         }
1011:         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_DIRICHLET_MARK;
1012:       }
1013:     }
1014:     ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);
1015:   }
1016:   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
1017:   if (graph->mirrors) {
1018:     for (i=0;i<graph->nvtxs;i++)
1019:       if (graph->mirrors[i])
1020:         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;

1022:     if (graph->xadj) {
1023:       PetscInt *new_xadj,*new_adjncy;
1024:       /* sort CSR graph */
1025:       for (i=0;i<graph->nvtxs;i++) {
1026:         PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);
1027:       }
1028:       /* adapt local CSR graph in case of local periodicity */
1029:       k = 0;
1030:       for (i=0;i<graph->nvtxs;i++)
1031:         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
1032:           k += graph->mirrors[graph->adjncy[j]];

1034:       PetscMalloc1(graph->nvtxs+1,&new_xadj);
1035:       PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);
1036:       new_xadj[0] = 0;
1037:       for (i=0;i<graph->nvtxs;i++) {
1038:         k = graph->xadj[i+1]-graph->xadj[i];
1039:         PetscArraycpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k);
1040:         new_xadj[i+1] = new_xadj[i]+k;
1041:         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
1042:           k = graph->mirrors[graph->adjncy[j]];
1043:           PetscArraycpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k);
1044:           new_xadj[i+1] += k;
1045:         }
1046:         k = new_xadj[i+1]-new_xadj[i];
1047:         PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);
1048:         new_xadj[i+1] = new_xadj[i]+k;
1049:       }
1050:       /* set new CSR into graph */
1051:       PetscFree(graph->xadj);
1052:       PetscFree(graph->adjncy);
1053:       graph->xadj = new_xadj;
1054:       graph->adjncy = new_adjncy;
1055:     }
1056:   }

1058:   /* mark special nodes (if any) -> each will become a single node equivalence class */
1059:   if (custom_primal_vertices) {
1060:     ISGetLocalSize(custom_primal_vertices,&is_size);
1061:     ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);
1062:     for (i=0,j=0;i<is_size;i++) {
1063:       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs  && graph->special_dof[is_indices[i]] != PCBDDCGRAPH_DIRICHLET_MARK) { /* out of bounds indices (if any) are skipped */
1064:         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_SPECIAL_MARK-j;
1065:         j++;
1066:       }
1067:     }
1068:     ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);
1069:   }

1071:   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
1072:   if (commsize > graph->commsizelimit) {
1073:     for (i=0;i<graph->nvtxs;i++) {
1074:       if (!graph->count[i]) {
1075:         PetscBTSet(graph->touched,i);
1076:         graph->subset[i] = 0;
1077:       }
1078:     }
1079:   }

1081:   /* init graph structure and compute default subsets */
1082:   nodes_touched = 0;
1083:   for (i=0;i<graph->nvtxs;i++) {
1084:     if (PetscBTLookup(graph->touched,i)) {
1085:       nodes_touched++;
1086:     }
1087:   }
1088:   i = 0;
1089:   graph->ncc = 0;
1090:   total_counts = 0;

1092:   /* allocated space for queues */
1093:   if (commsize == graph->commsizelimit) {
1094:     PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);
1095:   } else {
1096:     PetscInt nused = graph->nvtxs - nodes_touched;
1097:     PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);
1098:   }

1100:   while (nodes_touched<graph->nvtxs) {
1101:     /*  find first untouched node in local ordering */
1102:     while (PetscBTLookup(graph->touched,i)) i++;
1103:     PetscBTSet(graph->touched,i);
1104:     graph->subset[i] = graph->ncc+1;
1105:     graph->cptr[graph->ncc] = total_counts;
1106:     graph->queue[total_counts] = i;
1107:     total_counts++;
1108:     nodes_touched++;
1109:     /* now find all other nodes having the same set of sharing subdomains */
1110:     for (j=i+1;j<graph->nvtxs;j++) {
1111:       /* check for same number of sharing subdomains, dof number and same special mark */
1112:       if (!PetscBTLookup(graph->touched,j) && graph->count[i] == graph->count[j] && graph->which_dof[i] == graph->which_dof[j] && graph->special_dof[i] == graph->special_dof[j]) {
1113:         /* check for same set of sharing subdomains */
1114:         same_set = PETSC_TRUE;
1115:         for (k=0;k<graph->count[j];k++) {
1116:           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1117:             same_set = PETSC_FALSE;
1118:           }
1119:         }
1120:         /* I have found a friend of mine */
1121:         if (same_set) {
1122:           PetscBTSet(graph->touched,j);
1123:           graph->subset[j] = graph->ncc+1;
1124:           nodes_touched++;
1125:           graph->queue[total_counts] = j;
1126:           total_counts++;
1127:         }
1128:       }
1129:     }
1130:     graph->ncc++;
1131:   }
1132:   /* set default number of subsets (at this point no info on csr and/or local_subs has been taken into account, so n_subsets = ncc */
1133:   graph->n_subsets = graph->ncc;
1134:   PetscMalloc1(graph->n_subsets,&graph->subset_ncc);
1135:   for (i=0;i<graph->n_subsets;i++) {
1136:     graph->subset_ncc[i] = 1;
1137:   }
1138:   /* final pointer */
1139:   graph->cptr[graph->ncc] = total_counts;

1141:   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1142:   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1143:   PetscMalloc1(graph->ncc,&graph->subset_ref_node);
1144:   PetscMalloc1(graph->cptr[graph->ncc],&queue_global);
1145:   ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);
1146:   for (j=0;j<graph->ncc;j++) {
1147:     PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);
1148:     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1149:   }
1150:   PetscFree(queue_global);
1151:   graph->queue_sorted = PETSC_TRUE;

1153:   /* save information on subsets (needed when analyzing the connected components) */
1154:   if (graph->ncc) {
1155:     PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);
1156:     PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);
1157:     PetscArrayzero(graph->subset_idxs[0],graph->cptr[graph->ncc]);
1158:     for (j=1;j<graph->ncc;j++) {
1159:       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1160:       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1161:     }
1162:     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1163:     PetscArraycpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]);
1164:   }

1166:   /* renumber reference nodes */
1167:   ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);
1168:   ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);
1169:   ISDestroy(&subset_n);
1170:   ISRenumber(subset,NULL,NULL,&subset_n);
1171:   ISDestroy(&subset);
1172:   ISGetLocalSize(subset_n,&k);
1174:   ISGetIndices(subset_n,&is_indices);
1175:   PetscArraycpy(graph->subset_ref_node,is_indices,graph->ncc);
1176:   ISRestoreIndices(subset_n,&is_indices);
1177:   ISDestroy(&subset_n);

1179:   /* free workspace */
1180:   graph->setupcalled = PETSC_TRUE;
1181:   return 0;
1182: }

1184: PetscErrorCode PCBDDCGraphResetCoords(PCBDDCGraph graph)
1185: {
1186:   if (!graph) return 0;
1187:   PetscFree(graph->coords);
1188:   graph->cdim  = 0;
1189:   graph->cnloc = 0;
1190:   graph->cloc  = PETSC_FALSE;
1191:   return 0;
1192: }

1194: PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1195: {
1196:   if (!graph) return 0;
1197:   if (graph->freecsr) {
1198:     PetscFree(graph->xadj);
1199:     PetscFree(graph->adjncy);
1200:   } else {
1201:     graph->xadj = NULL;
1202:     graph->adjncy = NULL;
1203:   }
1204:   graph->freecsr = PETSC_FALSE;
1205:   graph->nvtxs_csr = 0;
1206:   return 0;
1207: }

1209: PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1210: {

1213:   if (!graph) return 0;
1214:   ISLocalToGlobalMappingDestroy(&graph->l2gmap);
1215:   PetscFree(graph->subset_ncc);
1216:   PetscFree(graph->subset_ref_node);
1217:   if (graph->nvtxs) {
1218:     PetscFree(graph->neighbours_set[0]);
1219:   }
1220:   PetscBTDestroy(&graph->touched);
1221:   PetscFree5(graph->count,
1222:                     graph->neighbours_set,
1223:                     graph->subset,
1224:                     graph->which_dof,
1225:                     graph->special_dof);
1226:   PetscFree2(graph->cptr,graph->queue);
1227:   if (graph->mirrors) {
1228:     PetscFree(graph->mirrors_set[0]);
1229:   }
1230:   PetscFree2(graph->mirrors,graph->mirrors_set);
1231:   if (graph->subset_idxs) {
1232:     PetscFree(graph->subset_idxs[0]);
1233:   }
1234:   PetscFree2(graph->subset_size,graph->subset_idxs);
1235:   ISDestroy(&graph->dirdofs);
1236:   ISDestroy(&graph->dirdofsB);
1237:   if (graph->n_local_subs) {
1238:     PetscFree(graph->local_subs);
1239:   }
1240:   graph->has_dirichlet       = PETSC_FALSE;
1241:   graph->twodimset           = PETSC_FALSE;
1242:   graph->twodim              = PETSC_FALSE;
1243:   graph->nvtxs               = 0;
1244:   graph->nvtxs_global        = 0;
1245:   graph->n_subsets           = 0;
1246:   graph->custom_minimal_size = 1;
1247:   graph->n_local_subs        = 0;
1248:   graph->maxcount            = PETSC_MAX_INT;
1249:   graph->setupcalled         = PETSC_FALSE;
1250:   return 0;
1251: }

1253: PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount)
1254: {
1255:   PetscInt       n;

1261:   /* raise an error if already allocated */
1263:   /* set number of vertices */
1264:   PetscObjectReference((PetscObject)l2gmap);
1265:   graph->l2gmap = l2gmap;
1266:   ISLocalToGlobalMappingGetSize(l2gmap,&n);
1267:   graph->nvtxs = n;
1268:   graph->nvtxs_global = N;
1269:   /* allocate used space */
1270:   PetscBTCreate(graph->nvtxs,&graph->touched);
1271:   PetscMalloc5(graph->nvtxs,&graph->count,graph->nvtxs,&graph->neighbours_set,graph->nvtxs,&graph->subset,graph->nvtxs,&graph->which_dof,graph->nvtxs,&graph->special_dof);
1272:   /* zeroes memory */
1273:   PetscArrayzero(graph->count,graph->nvtxs);
1274:   PetscArrayzero(graph->subset,graph->nvtxs);
1275:   /* use -1 as a default value for which_dof array */
1276:   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
1277:   PetscArrayzero(graph->special_dof,graph->nvtxs);
1278:   /* zeroes first pointer to neighbour set */
1279:   if (graph->nvtxs) {
1280:     graph->neighbours_set[0] = NULL;
1281:   }
1282:   /* zeroes workspace for values of ncc */
1283:   graph->subset_ncc = NULL;
1284:   graph->subset_ref_node = NULL;
1285:   /* maxcount for cc */
1286:   graph->maxcount = maxcount;
1287:   return 0;
1288: }

1290: PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1291: {
1292:   PCBDDCGraphResetCSR(*graph);
1293:   PCBDDCGraphResetCoords(*graph);
1294:   PCBDDCGraphReset(*graph);
1295:   PetscFree(*graph);
1296:   return 0;
1297: }

1299: PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1300: {
1301:   PCBDDCGraph    new_graph;

1303:   PetscNew(&new_graph);
1304:   new_graph->custom_minimal_size = 1;
1305:   new_graph->commsizelimit = 1;
1306:   *graph = new_graph;
1307:   return 0;
1308: }