Actual source code: Simplicializer.cxx

  1: #include <petscmesh.h>

  3: int debug = 0;

  5: #if 0

  9: PetscErrorCode ComputePreSievePartition(ALE::PreSieve* presieve, ALE::Point_set leaves, const char *name = NULL)
 10: {
 11:   MPI_Comm       comm = presieve->getComm();
 12:   PetscInt       numLeaves = leaves.size();
 13:   PetscMPIInt    rank, size;

 17:   ALE_LOG_EVENT_BEGIN
 18:   MPI_Comm_rank(comm, &rank);
 19:   MPI_Comm_size(comm, &size);
 20:   if (rank == 0) {
 21:     for(int p = 0; p < size; p++) {
 22:       ALE::Point partPoint(-1, p);
 23:       for(int l = (numLeaves/size)*p + PetscMin(numLeaves%size, p); l < (numLeaves/size)*(p+1) + PetscMin(numLeaves%size, p+1); l++) {
 24:         ALE::Point leaf(0, l);
 25:         ALE::Point_set cone = presieve->cone(leaf);
 26:         presieve->addCone(cone, partPoint);
 27:       }
 28:     }
 29:   } else {
 30:     ALE::Point partitionPoint(-1, rank);
 31:     presieve->addBasePoint(partitionPoint);
 32:   }
 33:   if (debug) {
 34:     ostringstream label1;
 35:     label1 << "Partition of presieve ";
 36:     if(name != NULL) {
 37:       label1 << "'" << name << "'";
 38:     }
 39:     label1 << "\n";
 40:     presieve->view(label1.str().c_str());
 41:   }
 42:   ALE_LOG_EVENT_END
 43:   return(0);
 44: }

 48: PetscErrorCode ComputeSievePartition(ALE::Sieve* sieve, const char *name = NULL)
 49: {
 50:   MPI_Comm       comm = sieve->getComm();
 51:   PetscInt       numLeaves = sieve->leaves().size();
 52:   PetscMPIInt    rank, size;

 56:   ALE_LOG_STAGE_BEGIN;
 57:   ALE_LOG_EVENT_BEGIN
 58:   MPI_Comm_rank(comm, &rank);
 59:   MPI_Comm_size(comm, &size);
 60:   if (rank == 0) {
 61:     for(int p = 0; p < size; p++) {
 62:       ALE::Point partPoint(-1, p);
 63:       for(int l = (numLeaves/size)*p + PetscMin(numLeaves%size, p); l < (numLeaves/size)*(p+1) + PetscMin(numLeaves%size, p+1); l++) {
 64:         sieve->addCone(sieve->closure(ALE::Point(0, l)), partPoint);
 65:       }
 66:     }
 67:   } else {
 68:     ALE::Point partitionPoint(-1, rank);
 69:     sieve->addBasePoint(partitionPoint);
 70:   }
 71:   if (debug) {
 72:     ostringstream label1;
 73:     label1 << "Partition of sieve ";
 74:     if(name != NULL) {
 75:       label1 << "'" << name << "'";
 76:     }
 77:     label1 << "\n";
 78:     sieve->view(label1.str().c_str());
 79:   }
 80:   ALE_LOG_EVENT_END
 81:   ALE_LOG_STAGE_END;
 82:   return(0);
 83: }

 87: PetscErrorCode PartitionPreSieve(ALE::PreSieve* presieve, const char *name = NULL, bool localize = 1, ALE::Obj<ALE::PreSieve> *pointTypes = NULL)
 88: {
 89:   MPI_Comm       comm = presieve->getComm();
 90:   PetscMPIInt    rank, size;

 94:   ALE_LOG_EVENT_BEGIN
 95:   MPI_Comm_rank(comm, &rank);
 96:   MPI_Comm_size(comm, &size);
 97:   // Cone complete to move the partitions to the other processors
 98:   ALE::Obj<ALE::Stack> completionStack = presieve->coneCompletion(ALE::PreSieve::completionTypePoint, ALE::PreSieve::footprintTypeCone, NULL);
 99:   ALE::Obj<ALE::PreSieve> completion = completionStack->top();
100:   if (debug) {
101:     ostringstream label1;
102:     label1 << "Completion";
103:     if(name != NULL) {
104:       label1 << " of '" << name << "'";
105:     }
106:     completion->view(label1.str().c_str());
107:   }
108:   // Create point type presieve
109:   if (pointTypes != NULL) {
110:     // Point types are as follows:
111:     //
112:     //   (rank, 0): Local point,  not shared with other processes
113:     //   (rank, 1): Leased point, owned but shared with other processes
114:     //     (otherRank, otherRank):    point leased to process otherRank
115:     //   (rank, 2): Rented point, not owned and shared with other processes
116:     //     (-otherRank-1, otherRank): point rented from process otherRank
117:     ALE::Obj<ALE::PreSieve> pTypes = ALE::PreSieve(comm);

119:     for(int p = 0; p < size; p++) {
120:       ALE::Point partitionPoint(-1, p);
121:       ALE::Point_set cone;

123:       if (presieve->baseContains(partitionPoint)) {
124:         cone = presieve->cone(partitionPoint);
125:       }

127:       if (p == rank) {
128:         ALE::Point point(rank, ALE::localPoint);
129:         pTypes->addCone(cone, point);

131:         if (completion->baseContains(partitionPoint)) {
132:           cone = completion->cone(partitionPoint);
133:         } else {
134:           cone.clear();
135:         }
136:         point = ALE::Point(rank, ALE::rentedPoint);
137:         pTypes->addCone(cone, point);
138:         for(ALE::Point_set::iterator e_itor = cone.begin(); e_itor != cone.end(); e_itor++) {
139:           ALE::Point e = *e_itor;
140:           ALE::Point f = *completionStack->support(e)->begin();

142:           point = ALE::Point(-f.index-1, f.index);
143:           pTypes->addCone(e, point);
144:         }
145:       } else {
146:         ALE::Point point;

148:         point = ALE::Point(rank, ALE::leasedPoint);
149:         pTypes->addCone(cone, point);
150:         point = ALE::Point(p, p);
151:         pTypes->addCone(cone, point);
152:       }
153:     }
154:     *pointTypes = pTypes;
155:     if (debug) {
156:       pTypes->view("Partition pointTypes");
157:     }
158:   }
159:   // Merge in the completion
160:   presieve->add(completion);
161:   // Move the cap to the base of the partition sieve
162:   ALE::Point partitionPoint(-1, rank);
163:   ALE::Point_set partition = presieve->cone(partitionPoint);
164:   for(ALE::Point_set::iterator p_itor = partition.begin(); p_itor != partition.end(); p_itor++) {
165:     ALE::Point p = *p_itor;
166:     presieve->addBasePoint(p);
167:   }
168:   if (debug) {
169:     ostringstream label2;
170:     if(name != NULL) {
171:       label2 << "Initial parallel state of '" << name << "'";
172:     } else {
173:       label2 << "Initial parallel presieve";
174:     }
175:     presieve->view(label2.str().c_str());
176:   }
177:   // Cone complete again to build the local topology
178:   completion = presieve->coneCompletion(ALE::PreSieve::completionTypePoint, ALE::PreSieve::footprintTypeCone, NULL)->top();
179:   if (debug) {
180:     ostringstream label3;
181:     if(name != NULL) {
182:       label3 << "Completion of '" << name << "'";
183:     } else {
184:       label3 << "Completion";
185:     }
186:     completion->view(label3.str().c_str());
187:   }
188:   presieve->add(completion);
189:   if (debug) {
190:     ostringstream label4;
191:     if(name != NULL) {
192:       label4 << "Completed parallel version of '" << name << "'";
193:     } else {
194:       label4 << "Completed parallel presieve";
195:     }
196:     presieve->view(label4.str().c_str());
197:   }
198:   // Unless explicitly prohibited, restrict to the local partition
199:   if(localize) {
200:     presieve->restrictBase(partition);
201:     if (debug) {
202:       ostringstream label5;
203:       if(name != NULL) {
204:         label5 << "Localized parallel version of '" << name << "'";
205:       } else {
206:         label5 << "Localized parallel presieve";
207:       }
208:       presieve->view(label5.str().c_str());
209:     }
210:   }
211:   // Support complete to get the adjacency information
212:   ALE_LOG_EVENT_END
213:   return(0);
214: }

218: /* This is currently duplicated in ex33mesh.c */
219: inline void ExpandInterval(ALE::Point interval, PetscInt indices[], PetscInt *indx)
220: {
221:   for(int i = 0; i < interval.index; i++) {
222:     indices[(*indx)++] = interval.prefix + i;
223:   }
224: }

228: /* This is currently duplicated in ex33mesh.c */
229: PetscErrorCode ExpandIntervals(ALE::Obj<ALE::Point_array> intervals, PetscInt *indices)
230: {
231:   int k = 0;

234:   for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
235:     ExpandInterval(*i_itor, indices, &k);
236:   }
237:   return(0);
238: }

242: PetscErrorCode MeshComputeGlobalScatter(ALE::Obj<ALE::IndexBundle> bundle, VecScatter *injection)
243: {
244:   VecScatter               sc;
245:   ALE::Obj<ALE::PreSieve>  globalIndices;
246:   ALE::Obj<ALE::PreSieve>  localIndices;
247:   ALE::Obj<ALE::Point_set> points;
248:   Vec                      l, g;
249:   IS                       localIS, globalIS;
250:   PetscInt                *localIdx, *globalIdx;
251:   PetscInt                 localSize, remoteSize, lcntr = 0, gcntr = 0;
252:   PetscErrorCode           ierr;

255:   localSize = bundle->getLocalSize();
256:   remoteSize = bundle->getRemoteSize();
257:   VecCreateSeq(PETSC_COMM_SELF, localSize+remoteSize, &l);
258:   VecCreateMPI(PETSC_COMM_WORLD, localSize, PETSC_DETERMINE, &g);
259:   PetscMalloc((localSize+remoteSize) * sizeof(PetscInt), &localIdx);
260:   PetscMalloc((localSize+remoteSize) * sizeof(PetscInt), &globalIdx);
261:   localIndices  = bundle->getLocalIndices();
262:   globalIndices = bundle->getGlobalIndices();
263:   points = globalIndices->base();
264:   for(ALE::Point_set::iterator p_itor = points->begin(); p_itor != points->end(); p_itor++) {
265:     ALE::Point p = *p_itor;
266:     ALE::Point_set lCone = localIndices->cone(p);
267:     ALE::Point_set gCone = globalIndices->cone(p);

269:     if (lCone.size()) {
270:       ExpandInterval(*(lCone.begin()), localIdx, &lcntr);
271:     }
272:     if (gCone.size()) {
273:       ExpandInterval(*(gCone.begin()), globalIdx, &gcntr);
274:     }
275:     if (lcntr != gcntr) {
276:       SETERRQ2(PETSC_ERR_PLIB, "Inconsistent numbering, %d != %d", lcntr, gcntr);
277:     }
278:   }
279:   if (lcntr != localSize+remoteSize) {
280:     SETERRQ2(PETSC_ERR_PLIB, "Inconsistent numbering, %d should be %d", lcntr, localSize+remoteSize);
281:   }
282:   ISCreateGeneral(PETSC_COMM_SELF, localSize+remoteSize, localIdx, &localIS);
283:   ISCreateGeneral(PETSC_COMM_SELF, localSize+remoteSize, globalIdx, &globalIS);
284:   PetscFree(localIdx);
285:   PetscFree(globalIdx);
286:   VecScatterCreate(l, localIS, g, globalIS, &sc);
287:   ISDestroy(localIS);
288:   ISDestroy(globalIS);
289:   *injection = sc;
290:   return(0);
291: }

295: PetscErrorCode setFiberValues(Vec b, ALE::Point e, ALE::IndexBundle* bundle, PetscScalar array[], InsertMode mode)
296: {
297:   ALE_LOG_EVENT_BEGIN
298:   ALE::Point_set   ee(e), empty;
299:   ALE::Obj<ALE::Point_set> intervals = bundle->getFiberIndices(ee, empty)->cap();
300:   static PetscInt  indicesSize = 0;
301:   static PetscInt *indices = NULL;
302:   PetscInt         numIndices = 0, i = 0;
303:   PetscErrorCode   ierr;

306:   for(ALE::Point_set::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
307:     numIndices += (*i_itor).index;
308:   }
309:   if (indicesSize && (indicesSize != numIndices)) {
310:     PetscFree(indices);
311:     indices = NULL;
312:   }
313:   if (!indices) {
314:     indicesSize = numIndices;
315:     PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
316:   }
317:   if (debug) {
318:     for(ALE::Point_set::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
319:       printf("indices (%d, %d)\n", (*i_itor).prefix, (*i_itor).index);
320:     }
321:   }
322:   ExpandInterval(*intervals->begin(), indices, &i);
323:   if (debug) {
324:     for(int i = 0; i < numIndices; i++) {
325:       printf("indices[%d] = %d\n", i, indices[i]);
326:     }
327:   }
328:   VecSetValues(b, numIndices, indices, array, mode);
329:   ALE_LOG_EVENT_END
330:   return(0);
331: }

335: PetscErrorCode setClosureValues(Vec b, ALE::Point e, ALE::IndexBundle* bundle, ALE::PreSieve* orientation, PetscScalar array[], InsertMode mode)
336: {
337:   ALE::Point_set   empty;
338:   ALE::Obj<ALE::Point_array> intervals = bundle->getClosureIndices(orientation->cone(e), empty);
339:   //ALE::Obj<ALE::Point_array> intervals = bundle->getOverlapOrderedIndices(orientation->cone(e), empty);
340:   static PetscInt  indicesSize = 0;
341:   static PetscInt *indices = NULL;
342:   PetscInt         numIndices = 0;
343:   PetscErrorCode   ierr;

346:   for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
347:     numIndices += (*i_itor).index;
348:   }
349:   if (indicesSize && (indicesSize != numIndices)) {
350:     PetscFree(indices);
351:     indices = NULL;
352:   }
353:   if (!indices) {
354:     indicesSize = numIndices;
355:     PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
356:   }
357:   for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
358:     printf("indices (%d, %d)\n", (*i_itor).prefix, (*i_itor).index);
359:   }
360:   ExpandIntervals(intervals, indices);
361:   for(int i = 0; i < numIndices; i++) {
362:     printf("indices[%d] = %d\n", i, indices[i]);
363:   }
364:   VecSetValues(b, numIndices, indices, array, mode);
365:   return(0);
366: }

370: PetscErrorCode restrictField(ALE::Obj<ALE::IndexBundle> bundle, ALE::Obj<ALE::PreSieve> orientation, PetscScalar *array, ALE::Point e, PetscScalar *values[])
371: {
372:   ALE::Obj<ALE::Point_array> intervals = bundle->getLocalOrderedClosureIndices(orientation->cone(e));
373:   /* This should be done by memory pooling by array size (we have a simple form below) */
374:   static PetscScalar *vals;
375:   static PetscInt     numValues = 0;
376:   static PetscInt    *indices = NULL;
377:   PetscInt            numIndices = 0;
378:   PetscErrorCode      ierr;

381:   for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
382:     numIndices += (*i_itor).index;
383:   }
384:   if (numValues && (numValues != numIndices)) {
385:     PetscFree(indices);
386:     indices = NULL;
387:     PetscFree(vals);
388:     vals = NULL;
389:   }
390:   if (!indices) {
391:     numValues = numIndices;
392:     PetscMalloc(numValues * sizeof(PetscInt), &indices);
393:     PetscMalloc(numValues * sizeof(PetscScalar), &vals);
394:   }
395:   if (debug) {
396:     for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
397:       printf("[%d]interval (%d, %d)\n", bundle->getCommRank(), (*i_itor).prefix, (*i_itor).index);
398:     }
399:   }
400:   ExpandIntervals(intervals, indices);
401:   for(int i = 0; i < numIndices; i++) {
402:     if (debug) {printf("[%d]indices[%d] = %d  val: %g\n", bundle->getCommRank(), i, indices[i], array[indices[i]]);}
403:     vals[i] = array[indices[i]];
404:   }
405:   *values = vals;
406:   return(0);
407: }

411: PetscErrorCode assembleField(ALE::Obj<ALE::IndexBundle> bundle, ALE::Obj<ALE::PreSieve> orientation, Vec b, ALE::Point e, PetscScalar array[], InsertMode mode)
412: {
413:   ALE::Obj<ALE::Point_array> intervals = bundle->getGlobalOrderedClosureIndices(orientation->cone(e));
414:   static PetscInt  indicesSize = 0;
415:   static PetscInt *indices = NULL;
416:   PetscInt         numIndices = 0;
417:   PetscErrorCode   ierr;

420:   for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
421:     numIndices += (*i_itor).index;
422:   }
423:   if (indicesSize && (indicesSize != numIndices)) {
424:     PetscFree(indices);
425:     indices = NULL;
426:   }
427:   if (!indices) {
428:     indicesSize = numIndices;
429:     PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
430:   }
431:   if (debug) {
432:     for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
433:       printf("[%d]Element (%d, %d) interval (%d, %d)\n", bundle->getCommRank(), e.prefix, e.index, (*i_itor).prefix, (*i_itor).index);
434:     }
435:   }
436:   ExpandIntervals(intervals, indices);
437:   if (debug) {
438:     for(int i = 0; i < numIndices; i++) {
439:       printf("[%d]indices[%d] = %d with value %g\n", bundle->getCommRank(), i, indices[i], array[i]);
440:     }
441:   }
442:   VecSetValues(b, numIndices, indices, array, mode);
443:   return(0);
444: }

448: PetscErrorCode assembleOperator(ALE::Obj<ALE::IndexBundle> bundle, ALE::Obj<ALE::PreSieve> orientation, Mat A, ALE::Point e, PetscScalar array[], InsertMode mode)
449: {
450:   ALE::Obj<ALE::Point_array> intervals = bundle->getGlobalOrderedClosureIndices(orientation->cone(e));
451:   static PetscInt  indicesSize = 0;
452:   static PetscInt *indices = NULL;
453:   PetscInt         numIndices = 0;
454:   PetscErrorCode   ierr;

457:   for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
458:     numIndices += (*i_itor).index;
459:   }
460:   if (indicesSize && (indicesSize != numIndices)) {
461:     PetscFree(indices);
462:     indices = NULL;
463:   }
464:   if (!indices) {
465:     indicesSize = numIndices;
466:     PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
467:   }
468:   if (debug) {
469:     for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
470:       printf("[%d]interval (%d, %d)\n", bundle->getCommRank(), (*i_itor).prefix, (*i_itor).index);
471:     }
472:   }
473:   ExpandIntervals(intervals, indices);
474:   if (debug) {
475:     for(int i = 0; i < numIndices; i++) {
476:       printf("[%d]indices[%d] = %d\n", bundle->getCommRank(), i, indices[i]);
477:     }
478:   }
479:   MatSetValues(A, numIndices, indices, numIndices, indices, array, mode);
480:   return(0);
481: }


488: PetscErrorCode createParallelCoordinates(Mesh mesh, int dim, ALE::Obj<ALE::PreSieve> partitionTypes)
489: {
490:   ALE::Obj<ALE::Sieve>       topology;
491:   ALE::Obj<ALE::IndexBundle> coordBundle;
492:   ALE::Obj<ALE::IndexBundle> serialCoordBundle;
493:   Vec                        coordinates, oldCoordinates, locCoordinates;
494:   VecScatter                 coordScatter;
495:   PetscErrorCode             ierr;

498:   ALE_LOG_EVENT_BEGIN
499:   MeshGetTopology(mesh, &topology);
500:   MeshGetCoordinateBundle(mesh, &serialCoordBundle);
501:   MeshGetCoordinates(mesh, &oldCoordinates);
502:   /* Create coordinate bundle */
503:   coordBundle = ALE::IndexBundle(topology);
504:   if (debug) {
505:     coordBundle->setVerbosity(11);
506:   }
507:   coordBundle->setFiberDimensionByDepth(0, dim);
508:   coordBundle->computeOverlapIndices();
509:   coordBundle->computeGlobalIndices();
510:   coordBundle->getLock();  // lock the bundle so that the overlap indices do not change
511:   MeshCreateVector(mesh, coordBundle, debug, &coordinates);
512:   VecSetBlockSize(coordinates, dim);
513:   VecGhostGetLocalForm(coordinates, &locCoordinates);
514:   VecSetBlockSize(locCoordinates, dim);
515:   MeshSetCoordinateBundle(mesh, coordBundle);
516:   /* Setup mapping to partitioned storage */
517:   MeshCreateMapping(mesh, serialCoordBundle, partitionTypes, coordBundle, &coordScatter);
518:   VecScatterBegin(oldCoordinates, coordinates, INSERT_VALUES, SCATTER_FORWARD, coordScatter);
519:   VecScatterEnd(oldCoordinates, coordinates, INSERT_VALUES, SCATTER_FORWARD, coordScatter);
520:   MeshSetCoordinates(mesh, coordinates);
521:   VecScatterDestroy(coordScatter);
522:   VecDestroy(oldCoordinates);
523:   if (debug) {
524:     MPI_Comm comm;
525:     PetscObjectGetComm((PetscObject) mesh, &comm);
526:     PetscPrintf(comm, "Parallel Coordinates\n===========================\n");
527:     VecView(coordinates, PETSC_VIEWER_STDOUT_WORLD);
528:   }
529:   /* Communicate ghosted coordinates */
530:   VecGhostUpdateBegin(coordinates, INSERT_VALUES, SCATTER_FORWARD);
531:   VecGhostUpdateEnd(coordinates, INSERT_VALUES, SCATTER_FORWARD);
532:   ALE_LOG_EVENT_END
533:   return(0);
534: }

538: PetscErrorCode MeshCreateMapping(Mesh mesh, ALE::Obj<ALE::IndexBundle> sourceBundle, ALE::Obj<ALE::PreSieve> pointTypes, ALE::Obj<ALE::IndexBundle> targetBundle, VecScatter *scatter)
539: {
540:   ALE::Obj<ALE::Stack>     mappingStack;
541:   ALE::Obj<ALE::PreSieve>  sourceIndices, targetIndices;
542:   ALE::Obj<ALE::Point_set> base;
543:   Vec                      sourceVec, targetVec;
544:   IS                       fromIS, toIS;
545:   PetscInt                *fromIndices, *toIndices;
546:   PetscInt                 fromIdx = 0, toIdx = 0;
547:   PetscInt                 locSourceSize = sourceBundle->getLocalSize();
548:   PetscErrorCode           ierr;

551:   mappingStack  = sourceBundle->computeMappingIndices(pointTypes, targetBundle);
552:   sourceIndices = mappingStack->top();
553:   targetIndices = mappingStack->bottom();
554:   base = sourceIndices->base();
555: #if 1
556:   locSourceSize = 0;
557:   for(ALE::Point_set::iterator e_itor = base->begin(); e_itor != base->end(); e_itor++) {
558:     ALE::Obj<ALE::Point_set> sourceCone = sourceIndices->cone(*e_itor);
559:     ALE::Obj<ALE::Point_set> targetCone = targetIndices->cone(*e_itor);

561:     if (sourceCone->size() && targetCone->size()) {
562:       if (sourceCone->begin()->index != targetCone->begin()->index) {
563:         SETERRQ2(PETSC_ERR_PLIB, "Mismatch in index sizes %d and %d", sourceCone->begin()->index, targetCone->begin()->index);
564:       }
565:       locSourceSize += sourceCone->begin()->index;
566:     }
567:   }
568: #endif
569:   PetscMalloc(locSourceSize * sizeof(PetscInt), &fromIndices);
570:   PetscMalloc(locSourceSize * sizeof(PetscInt), &toIndices);
571:   for(ALE::Point_set::iterator e_itor = base->begin(); e_itor != base->end(); e_itor++) {
572:     ALE::Obj<ALE::Point_set> sourceCone = sourceIndices->cone(*e_itor);
573:     ALE::Obj<ALE::Point_set> targetCone = targetIndices->cone(*e_itor);

575:     if (sourceCone->size() && targetCone->size()) {
576:       ExpandInterval(*sourceCone->begin(), fromIndices, &fromIdx);
577:       ExpandInterval(*targetCone->begin(), toIndices,   &toIdx);
578:     }
579:   }
580: #if 0
581:   if ((fromIdx != locSourceSize) || (toIdx != locSourceSize)) {
582:     SETERRQ3(PETSC_ERR_PLIB, "Invalid index sizes %d, %d should be %d", fromIdx, toIdx, locSourceSize);
583:   }
584: #endif
585:   ISCreateGeneral(PETSC_COMM_SELF, locSourceSize, fromIndices, &fromIS);
586:   ISCreateGeneral(PETSC_COMM_SELF, locSourceSize, toIndices,   &toIS);
587:   PetscFree(fromIndices);
588:   PetscFree(toIndices);
589:   MeshCreateVector(mesh, sourceBundle, debug, &sourceVec);
590:   MeshCreateVector(mesh, targetBundle, debug, &targetVec);
591:   VecScatterCreate(sourceVec, fromIS, targetVec, toIS, scatter);
592:   VecDestroy(sourceVec);
593:   VecDestroy(targetVec);
594:   ISDestroy(fromIS);
595:   ISDestroy(toIS);
596:   return(0);
597: }

601: /*@
602:   MeshDistribute - 
603: */
604: PetscErrorCode MeshDistribute(Mesh mesh)
605: {
606:   ALE::Obj<ALE::Sieve>       topology;
607:   ALE::Obj<ALE::PreSieve>    orientation;
608:   ALE::Obj<ALE::IndexBundle> elementBundle;
609:   ALE::Obj<ALE::PreSieve>    partitionTypes;
610:   MPI_Comm                   comm;
611:   PetscMPIInt                rank;
612:   PetscInt                   dim;
613:   PetscErrorCode             ierr;

617:   ALE_LOG_EVENT_BEGIN
618:   PetscObjectGetComm((PetscObject) mesh, &comm);
619:   MPI_Comm_rank(comm, &rank);
620:   MeshGetTopology(mesh, &topology);
621:   MeshGetOrientation(mesh, &orientation);
622:   dim = topology->diameter();
623:   topology->setStratificationPolicy(ALE::Sieve::stratificationPolicyOnLocking);
624:   /* Partition the topology and orientation */
625:   ComputePreSievePartition(orientation, topology->leaves(), "Orientation");
626:   PartitionPreSieve(orientation, "Orientation", 1);
627:   ComputeSievePartition(topology, "Topology");
628:   PartitionPreSieve(topology, "Topology", 1, &partitionTypes);
629:   topology->getLock();
630:   topology->releaseLock();
631:   topology->setStratificationPolicy(ALE::Sieve::stratificationPolicyOnMutation);
632:   /* Add the trivial vertex orientation */
633:   ALE::Obj<ALE::Point_set> roots = topology->depthStratum(0);
634:   for(ALE::Point_set::iterator vertex_itor = roots->begin(); vertex_itor != roots->end(); vertex_itor++) {
635:     ALE::Point v = *vertex_itor;
636:     orientation->addCone(v, v);
637:   }
638:   /* Create element bundle */
639:   elementBundle = ALE::IndexBundle(topology);
640:   if (debug) {
641:     elementBundle->setVerbosity(11);
642:   }
643:   elementBundle->setFiberDimensionByHeight(0, 1);
644:   elementBundle->computeOverlapIndices();
645:   elementBundle->computeGlobalIndices();
646:   elementBundle->getLock();  // lock the bundle so that the overlap indices do not change
647:   MeshSetElementBundle(mesh, elementBundle);
648:   createParallelCoordinates(mesh, dim, partitionTypes);
649:   ALE_LOG_EVENT_END
650:   return(0);
651: }

655: /*@
656:   MeshUnify - 
657: */
658: PetscErrorCode MeshUnify(Mesh mesh, Mesh *serialMesh)
659: {
660:   ALE::Obj<ALE::Sieve>       topology;
661:   ALE::Obj<ALE::PreSieve>    orientation;
662:   ALE::Obj<ALE::IndexBundle> coordBundle;
663:   ALE::Obj<ALE::PreSieve>    partitionTypes;
664:   ALE::Obj<ALE::Sieve>       boundary;
665:   Vec               serialCoordinates, coordinates;
666:   MPI_Comm          comm;
667:   PetscMPIInt       rank;
668:   PetscInt          dim;
669:   PetscErrorCode    ierr;

673:   PetscObjectGetComm((PetscObject) mesh, &comm);
674:   MPI_Comm_rank(comm, &rank);
675:   MeshGetCoordinateBundle(mesh, &coordBundle);
676:   MeshGetCoordinates(mesh, &coordinates);
677:   /* Unify the topology and orientation */
678:   MeshGetTopology(mesh, &topology);
679:   MeshGetOrientation(mesh, &orientation);
680:   ALE::Obj<ALE::Sieve>     serialTopology = ALE::Sieve(comm);
681:   ALE::Obj<ALE::PreSieve>  serialOrientation = ALE::PreSieve(comm);
682:   ALE::Obj<ALE::Point_set> base = topology->base();
683:   ALE::Obj<ALE::Point_set> orientationBase = orientation->base();
684:   ALE::Point               partition(-1, 0);

686:   MeshCreate(comm, serialMesh);
687:   serialOrientation->addCone(orientation->space(), partition);
688:   for(ALE::Point_set::iterator b_itor = orientationBase->begin(); b_itor != orientationBase->end(); b_itor++) {
689:     serialOrientation->addCone(orientation->cone(*b_itor), *b_itor);
690:   }
691:   PartitionPreSieve(serialOrientation, "Serial Orientation", 1);
692:   serialTopology->addCone(topology->space(), partition);
693:   for(ALE::Point_set::iterator b_itor = base->begin(); b_itor != base->end(); b_itor++) {
694:     serialTopology->addCone(topology->cone(*b_itor), *b_itor);
695:   }
696:   PartitionPreSieve(serialTopology, "Serial Topology", 1, &partitionTypes);
697:   MeshSetTopology(*serialMesh, serialTopology);
698:   MeshSetOrientation(*serialMesh, serialOrientation);
699:   /* Unify boundary */
700:   MeshGetBoundary(mesh, &boundary);
701:   ALE::Obj<ALE::Sieve> serialBoundary = ALE::Sieve(comm);
702:   ALE::Obj<ALE::Point_set> boundaryBase = boundary->base();

704:   serialBoundary->addCone(boundary->space(), partition);
705:   for(ALE::Point_set::iterator b_itor = boundaryBase->begin(); b_itor != boundaryBase->end(); b_itor++) {
706:     serialBoundary->addCone(boundary->cone(*b_itor), *b_itor);
707:   }
708:   PartitionPreSieve(serialBoundary, "Serial Boundary", 1);
709:   MeshSetBoundary(*serialMesh, serialBoundary);
710:   /* Create vertex bundle */
711:   ALE::Obj<ALE::IndexBundle> serialVertexBundle = ALE::IndexBundle(serialTopology);
712:   serialVertexBundle->setFiberDimensionByDepth(0, 1);
713:   serialVertexBundle->computeOverlapIndices();
714:   serialVertexBundle->computeGlobalIndices();
715:   serialVertexBundle->getLock();  // lock the bundle so that the overlap indices do not change
716:   MeshSetVertexBundle(*serialMesh, serialVertexBundle);
717:   /* Create element bundle */
718:   ALE::Obj<ALE::IndexBundle> serialElementBundle = ALE::IndexBundle(serialTopology);
719:   serialElementBundle->setFiberDimensionByHeight(0, 1);
720:   serialElementBundle->computeOverlapIndices();
721:   serialElementBundle->computeGlobalIndices();
722:   serialElementBundle->getLock();  // lock the bundle so that the overlap indices do not change
723:   MeshSetElementBundle(*serialMesh, serialElementBundle);
724:   /* Create coordinate bundle and storage */
725:   ALE::Obj<ALE::IndexBundle> serialCoordBundle = ALE::IndexBundle(serialTopology);
726:   dim = topology->diameter();
727:   serialCoordBundle->setFiberDimensionByDepth(0, dim);
728:   serialCoordBundle->computeOverlapIndices();
729:   serialCoordBundle->computeGlobalIndices();
730:   serialCoordBundle->getLock();  // lock the bundle so that the overlap indices do not change
731:   MeshCreateVector(*serialMesh, serialCoordBundle, debug, &serialCoordinates);
732:   VecSetBlockSize(serialCoordinates, dim);
733:   MeshSetCoordinateBundle(*serialMesh, serialCoordBundle);
734:   /* Setup mapping to unified storage */
735:   VecScatter coordScatter;

737:   MeshCreateMapping(mesh, coordBundle, partitionTypes, serialCoordBundle, &coordScatter);
738:   VecScatterBegin(coordinates, serialCoordinates, INSERT_VALUES, SCATTER_FORWARD, coordScatter);
739:   VecScatterEnd(coordinates, serialCoordinates, INSERT_VALUES, SCATTER_FORWARD, coordScatter);
740:   MeshSetCoordinates(*serialMesh, serialCoordinates);
741:   VecScatterDestroy(coordScatter);
742:   return(0);
743: }

745: ALE::Obj<ALE::Point_set> getLocal(MPI_Comm comm, ALE::Obj<ALE::Stack> spaceFootprint, ALE::Obj<ALE::Point_set> points)
746: {
747:   ALE::Obj<ALE::Point_set> localPoints(new ALE::Point_set);
748:   ALE::Point     proc(0, spaceFootprint->getCommRank());

750:   for(ALE::Point_set::iterator p_itor = points->begin(); p_itor != points->end(); p_itor++) {
751:     if (*spaceFootprint->cone(*p_itor)->begin() != proc) continue;
752:     localPoints->insert(*p_itor);
753:   }
754:   return localPoints;
755: }

757: PetscErrorCode MeshComputeOverlap(Mesh mesh)
758: {
759:   ALE::Obj<ALE::Sieve> topology;
760:   ALE::Obj<ALE::Stack> spaceFootprint;
761:   PetscErrorCode       ierr;

764:   MeshGetSpaceFootprint(mesh, &spaceFootprint);
765:   if (!spaceFootprint) {
766:     MeshGetTopology(mesh, &topology);
767:     spaceFootprint = topology->spaceFootprint(ALE::PreSieve::completionTypePoint, ALE::PreSieve::footprintTypeSupport, NULL);
768:     MeshSetSpaceFootprint(mesh, spaceFootprint);
769:   }
770:   return(0);
771: }

773: PetscErrorCode MeshGetDimLocalSize(Mesh mesh, int dim, PetscInt *size)
774: {
775:   MPI_Comm             comm;
776:   ALE::Obj<ALE::Sieve> topology;
777:   ALE::Obj<ALE::Stack> spaceFootprint;
778:   PetscErrorCode       ierr;

783:   PetscObjectGetComm((PetscObject) mesh, &comm);
784:   MeshComputeOverlap(mesh);
785:   MeshGetTopology(mesh, &topology);
786:   MeshGetSpaceFootprint(mesh, &spaceFootprint);
787:   *size = getLocal(comm, *spaceFootprint, topology->depthStratum(dim))->size();
788:   return(0);
789: }

791: PetscErrorCode MeshGetDimLocalRanges(Mesh mesh, int dim, PetscInt starts[])
792: {
793:   MPI_Comm       comm;
794:   PetscInt       localSize;

800:   PetscObjectGetComm((PetscObject) mesh, &comm);
801:   MeshGetDimLocalSize(mesh, dim, &localSize);
802:   MPI_Allgather(&localSize, 1, MPI_INT, starts, 1, MPI_INT, comm);
803:   return(0);
804: }

806: PetscErrorCode MeshGetDimGlobalSize(Mesh mesh, int dim, PetscInt *size)
807: {
808:   MPI_Comm       comm;
809:   PetscInt       localSize, globalSize;

815:   PetscObjectGetComm((PetscObject) mesh, &comm);
816:   MeshGetDimLocalSize(mesh, dim, &localSize);
817:   MPI_Allreduce(&localSize, &globalSize, 1, MPI_INT, MPI_SUM, comm);
818:   *size = globalSize;
819:   return(0);
820: }

822: #endif