Actual source code: Mesh.hh

  1: #ifndef included_ALE_Mesh_hh
  2: #define included_ALE_Mesh_hh

  4: #ifndef  included_ALE_Numbering_hh
  5: #include <Numbering.hh>
  6: #endif

  8: #ifndef  included_ALE_Field_hh
  9: #include <Field.hh>
 10: #endif

 12: #ifndef  included_ALE_SieveBuilder_hh
 13: #include <SieveBuilder.hh>
 14: #endif

 16: namespace ALE {
 17:   template<typename Sieve_,
 18:            typename RealSection_  = Section<typename Sieve_::point_type, double>,
 19:            typename IntSection_   = Section<typename Sieve_::point_type, int>,
 20:            typename ArrowSection_ = UniformSection<MinimalArrow<typename Sieve_::point_type, typename Sieve_::point_type>, int> >
 21:   class Bundle : public ALE::ParallelObject {
 22:   public:
 23:     typedef Sieve_                                                    sieve_type;
 24:     typedef RealSection_                                              real_section_type;
 25:     typedef IntSection_                                               int_section_type;
 26:     typedef ArrowSection_                                             arrow_section_type;
 27:     typedef Bundle<Sieve_,RealSection_,IntSection_,ArrowSection_>     this_type;
 28:     typedef typename sieve_type::point_type                           point_type;
 29:     typedef typename ALE::Sifter<int, point_type, int>                label_type;
 30:     typedef typename std::map<const std::string, Obj<label_type> >    labels_type;
 31:     typedef typename label_type::supportSequence                      label_sequence;
 32:     typedef std::map<std::string, Obj<arrow_section_type> >           arrow_sections_type;
 33:     typedef std::map<std::string, Obj<real_section_type> >            real_sections_type;
 34:     typedef std::map<std::string, Obj<int_section_type> >             int_sections_type;
 35:     typedef ALE::Point                                                index_type;
 36:     typedef std::pair<index_type, int>                                oIndex_type;
 37:     typedef std::vector<oIndex_type>                                  oIndexArray;
 38:     typedef std::pair<int *, int>                                     indices_type;
 39:     typedef NumberingFactory<this_type>                               NumberingFactory;
 40:     typedef typename NumberingFactory::numbering_type                 numbering_type;
 41:     typedef typename NumberingFactory::order_type                     order_type;
 42:     typedef typename ALE::Sifter<int,point_type,point_type>           send_overlap_type;
 43:     typedef typename ALE::Sifter<point_type,int,point_type>           recv_overlap_type;
 44:     typedef typename ALE::SieveAlg<this_type>                         sieve_alg_type;
 45:     typedef typename sieve_alg_type::coneArray                        coneArray;
 46:     typedef typename sieve_alg_type::orientedConeArray                oConeArray;
 47:     typedef typename sieve_alg_type::supportArray                     supportArray;
 48:   protected:
 49:     Obj<sieve_type>       _sieve;
 50:     labels_type           _labels;
 51:     int                   _maxHeight;
 52:     int                   _maxDepth;
 53:     arrow_sections_type   _arrowSections;
 54:     real_sections_type    _realSections;
 55:     int_sections_type     _intSections;
 56:     Obj<oIndexArray>      _indexArray;
 57:     Obj<NumberingFactory> _factory;
 58:     bool                   _calculatedOverlap;
 59:     Obj<send_overlap_type> _sendOverlap;
 60:     Obj<recv_overlap_type> _recvOverlap;
 61:     Obj<send_overlap_type> _distSendOverlap;
 62:     Obj<recv_overlap_type> _distRecvOverlap;
 63:     // Work space
 64:     Obj<std::set<point_type> > _modifiedPoints;
 65:   public:
 66:     Bundle(MPI_Comm comm, int debug = 0) : ALE::ParallelObject(comm, debug), _maxHeight(-1), _maxDepth(-1) {
 67:       this->_indexArray        = new oIndexArray();
 68:       this->_modifiedPoints    = new std::set<point_type>();
 69:       this->_factory           = NumberingFactory::singleton(this->comm(), this->debug());
 70:       this->_calculatedOverlap = false;
 71:       this->_sendOverlap       = new send_overlap_type(comm, debug);
 72:       this->_recvOverlap       = new recv_overlap_type(comm, debug);
 73:     };
 74:     Bundle(const Obj<sieve_type>& sieve) : ALE::ParallelObject(sieve->comm(), sieve->debug()), _sieve(sieve), _maxHeight(-1), _maxDepth(-1) {
 75:       this->_indexArray        = new oIndexArray();
 76:       this->_modifiedPoints    = new std::set<point_type>();
 77:       this->_factory           = NumberingFactory::singleton(this->comm(), this->debug());
 78:       this->_calculatedOverlap = false;
 79:       this->_sendOverlap       = new send_overlap_type(comm, debug);
 80:       this->_recvOverlap       = new recv_overlap_type(comm, debug);
 81:     };
 82:     virtual ~Bundle() {};
 83:   public: // Verifiers
 84:     bool hasLabel(const std::string& name) {
 85:       if (this->_labels.find(name) != this->_labels.end()) {
 86:         return true;
 87:       }
 88:       return false;
 89:     };
 90:     void checkLabel(const std::string& name) {
 91:       if (!this->hasLabel(name)) {
 92:         ostringstream msg;
 93:         msg << "Invalid label name: " << name << std::endl;
 94:         throw ALE::Exception(msg.str().c_str());
 95:       }
 96:     };
 97:   public: // Accessors
 98:     const Obj<sieve_type>& getSieve() const {return this->_sieve;};
 99:     void setSieve(const Obj<sieve_type>& sieve) {this->_sieve = sieve;};
100:     bool hasArrowSection(const std::string& name) const {
101:       return this->_arrowSections.find(name) != this->_arrowSections.end();
102:     };
103:     const Obj<arrow_section_type>& getArrowSection(const std::string& name) {
104:       if (!this->hasArrowSection(name)) {
105:         Obj<arrow_section_type> section = new arrow_section_type(this->comm(), this->debug());

107:         section->setName(name);
108:         if (this->_debug) {std::cout << "Creating new arrow section: " << name << std::endl;}
109:         this->_arrowSections[name] = section;
110:       }
111:       return this->_arrowSections[name];
112:     };
113:     void setArrowSection(const std::string& name, const Obj<arrow_section_type>& section) {
114:       this->_arrowSections[name] = section;
115:     };
116:     Obj<std::set<std::string> > getArrowSections() const {
117:       Obj<std::set<std::string> > names = std::set<std::string>();

119:       for(typename arrow_sections_type::const_iterator s_iter = this->_arrowSections.begin(); s_iter != this->_arrowSections.end(); ++s_iter) {
120:         names->insert(s_iter->first);
121:       }
122:       return names;
123:     };
124:     bool hasRealSection(const std::string& name) const {
125:       return this->_realSections.find(name) != this->_realSections.end();
126:     };
127:     const Obj<real_section_type>& getRealSection(const std::string& name) {
128:       if (!this->hasRealSection(name)) {
129:         Obj<real_section_type> section = new real_section_type(this->comm(), this->debug());

131:         section->setName(name);
132:         if (this->_debug) {std::cout << "Creating new real section: " << name << std::endl;}
133:         this->_realSections[name] = section;
134:       }
135:       return this->_realSections[name];
136:     };
137:     void setRealSection(const std::string& name, const Obj<real_section_type>& section) {
138:       this->_realSections[name] = section;
139:     };
140:     Obj<std::set<std::string> > getRealSections() const {
141:       Obj<std::set<std::string> > names = std::set<std::string>();

143:       for(typename real_sections_type::const_iterator s_iter = this->_realSections.begin(); s_iter != this->_realSections.end(); ++s_iter) {
144:         names->insert(s_iter->first);
145:       }
146:       return names;
147:     };
148:     bool hasIntSection(const std::string& name) const {
149:       return this->_intSections.find(name) != this->_intSections.end();
150:     };
151:     const Obj<int_section_type>& getIntSection(const std::string& name) {
152:       if (!this->hasIntSection(name)) {
153:         Obj<int_section_type> section = new int_section_type(this->comm(), this->debug());

155:         section->setName(name);
156:         if (this->_debug) {std::cout << "Creating new int section: " << name << std::endl;}
157:         this->_intSections[name] = section;
158:       }
159:       return this->_intSections[name];
160:     };
161:     void setIntSection(const std::string& name, const Obj<int_section_type>& section) {
162:       this->_intSections[name] = section;
163:     };
164:     Obj<std::set<std::string> > getIntSections() const {
165:       Obj<std::set<std::string> > names = std::set<std::string>();

167:       for(typename int_sections_type::const_iterator s_iter = this->_intSections.begin(); s_iter != this->_intSections.end(); ++s_iter) {
168:         names->insert(s_iter->first);
169:       }
170:       return names;
171:     };
172:     const Obj<NumberingFactory>& getFactory() const {return this->_factory;};
173:     const Obj<send_overlap_type>& getSendOverlap() const {return this->_sendOverlap;};
174:     void setSendOverlap(const Obj<send_overlap_type>& overlap) {this->_sendOverlap = overlap;};
175:     const Obj<recv_overlap_type>& getRecvOverlap() const {return this->_recvOverlap;};
176:     void setRecvOverlap(const Obj<recv_overlap_type>& overlap) {this->_recvOverlap = overlap;};
177:     const Obj<send_overlap_type>& getDistSendOverlap() const {return this->_distSendOverlap;};
178:     void setDistSendOverlap(const Obj<send_overlap_type>& overlap) {this->_distSendOverlap = overlap;};
179:     const Obj<recv_overlap_type>& getDistRecvOverlap() const {return this->_distRecvOverlap;};
180:     void setDistRecvOverlap(const Obj<recv_overlap_type>& overlap) {this->_distRecvOverlap = overlap;};
181:   public: // Labels
182:     int getValue (const Obj<label_type>& label, const point_type& point, const int defValue = 0) {
183:       const Obj<typename label_type::coneSequence>& cone = label->cone(point);

185:       if (cone->size() == 0) return defValue;
186:       return *cone->begin();
187:     };
188:     void setValue(const Obj<label_type>& label, const point_type& point, const int value) {
189:       label->setCone(value, point);
190:     };
191:     template<typename InputPoints>
192:     int getMaxValue (const Obj<label_type>& label, const Obj<InputPoints>& points, const int defValue = 0) {
193:       int maxValue = defValue;

195:       for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
196:         maxValue = std::max(maxValue, this->getValue(label, *p_iter, defValue));
197:       }
198:       return maxValue;
199:     };
200:     const Obj<label_type>& createLabel(const std::string& name) {
201:       this->_labels[name] = new label_type(this->comm(), this->debug());
202:       return this->_labels[name];
203:     };
204:     const Obj<label_type>& getLabel(const std::string& name) {
205:       this->checkLabel(name);
206:       return this->_labels[name];
207:     };
208:     const labels_type& getLabels() {
209:       return this->_labels;
210:     };
211:     virtual const Obj<label_sequence>& getLabelStratum(const std::string& name, int value) {
212:       this->checkLabel(name);
213:       return this->_labels[name]->support(value);
214:     };
215:   public: // Stratification
216:     template<class InputPoints>
217:     void computeHeight(const Obj<label_type>& height, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxHeight) {
218:       this->_modifiedPoints->clear();

220:       for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
221:         // Compute the max height of the points in the support of p, and add 1
222:         int h0 = this->getValue(height, *p_iter, -1);
223:         int h1 = this->getMaxValue(height, sieve->support(*p_iter), -1) + 1;

225:         if(h1 != h0) {
226:           this->setValue(height, *p_iter, h1);
227:           if (h1 > maxHeight) maxHeight = h1;
228:           this->_modifiedPoints->insert(*p_iter);
229:         }
230:       }
231:       // FIX: We would like to avoid the copy here with cone()
232:       if(this->_modifiedPoints->size() > 0) {
233:         this->computeHeight(height, sieve, sieve->cone(this->_modifiedPoints), maxHeight);
234:       }
235:     };
236:     void computeHeights() {
237:       const Obj<label_type>& label = this->createLabel(std::string("height"));

239:       this->_maxHeight = -1;
240:       this->computeHeight(label, this->_sieve, this->_sieve->leaves(), this->_maxHeight);
241:     };
242:     virtual int height() const {return this->_maxHeight;};
243:     virtual int height(const point_type& point) {
244:       return this->getValue(this->_labels["height"], point, -1);
245:     };
246:     virtual const Obj<label_sequence>& heightStratum(int height) {
247:       return this->getLabelStratum("height", height);
248:     };
249:     template<class InputPoints>
250:     void computeDepth(const Obj<label_type>& depth, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxDepth) {
251:       this->_modifiedPoints->clear();

253:       for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
254:         // Compute the max depth of the points in the cone of p, and add 1
255:         int d0 = this->getValue(depth, *p_iter, -1);
256:         int d1 = this->getMaxValue(depth, sieve->cone(*p_iter), -1) + 1;

258:         if(d1 != d0) {
259:           this->setValue(depth, *p_iter, d1);
260:           if (d1 > maxDepth) maxDepth = d1;
261:           this->_modifiedPoints->insert(*p_iter);
262:         }
263:       }
264:       // FIX: We would like to avoid the copy here with support()
265:       if(this->_modifiedPoints->size() > 0) {
266:         this->computeDepth(depth, sieve, sieve->support(this->_modifiedPoints), maxDepth);
267:       }
268:     };
269:     void computeDepths() {
270:       const Obj<label_type>& label = this->createLabel(std::string("depth"));

272:       this->_maxDepth = -1;
273:       this->computeDepth(label, this->_sieve, this->_sieve->roots(), this->_maxDepth);
274:     };
275:     virtual int depth() const {return this->_maxDepth;};
276:     virtual int depth(const point_type& point) {
277:       return this->getValue(this->_labels["depth"], point, -1);
278:     };
279:     virtual const Obj<label_sequence>& depthStratum(int depth) {
280:       return this->getLabelStratum("depth", depth);
281:     };
284:     virtual void stratify() {
285:       ALE_LOG_EVENT_BEGIN;
286:       this->computeHeights();
287:       this->computeDepths();
288:       ALE_LOG_EVENT_END;
289:     };
290:   public: // Size traversal
291:     template<typename Section_>
292:     int size(const Obj<Section_>& section, const point_type& p) {
293:       const typename Section_::chart_type&  chart   = section->getChart();
294:       const Obj<coneArray>                  closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), p);
295:       typename coneArray::iterator          end     = closure->end();
296:       int                                   size    = 0;

298:       for(typename coneArray::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
299:         if (chart.count(*c_iter)) {
300:           size += section->getConstrainedFiberDimension(*c_iter);
301:         }
302:       }
303:       return size;
304:     };
305:     template<typename Section_>
306:     int sizeWithBC(const Obj<Section_>& section, const point_type& p) {
307:       const typename Section_::chart_type&  chart   = section->getChart();
308:       const Obj<coneArray>                  closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), p);
309:       typename coneArray::iterator          end     = closure->end();
310:       int                                   size    = 0;

312:       for(typename coneArray::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
313:         if (chart.count(*c_iter)) {
314:           size += section->getFiberDimension(*c_iter);
315:         }
316:       }
317:       return size;
318:     };
319:   protected:
320:     int *getIndexArray(const int size) {
321:       static int *array   = NULL;
322:       static int  maxSize = 0;

324:       if (size > maxSize) {
325:         maxSize = size;
326:         if (array) delete [] array;
327:         array = new int[maxSize];
328:       };
329:       return array;
330:     };
331:   public: // Index traversal
332:     void expandInterval(const index_type& interval, PetscInt indices[], PetscInt *indx) {
333:       const int end = interval.prefix + interval.index;

335:       for(int i = interval.index; i < end; ++i) {
336:         indices[(*indx)++] = i;
337:       }
338:     };
339:     void expandInterval(const index_type& interval, const int orientation, PetscInt indices[], PetscInt *indx) {
340:       if (orientation >= 0) {
341:         for(int i = 0; i < interval.prefix; ++i) {
342:           indices[(*indx)++] = interval.index + i;
343:         }
344:       } else {
345:         for(int i = interval.prefix-1; i >= 0; --i) {
346:           indices[(*indx)++] = interval.index + i;
347:         }
348:       }
349:       for(int i = 0; i < -interval.prefix; ++i) {
350:         indices[(*indx)++] = -1;
351:       }
352:     };
353:     void expandIntervals(Obj<oIndexArray> intervals, PetscInt *indices) {
354:       int k = 0;

356:       for(typename oIndexArray::iterator i_iter = intervals->begin(); i_iter != intervals->end(); ++i_iter) {
357:         this->expandInterval(i_iter->first, i_iter->second, indices, &k);
358:       }
359:     }
360:     template<typename Section_>
361:     const indices_type getIndices(const Obj<Section_>& section, const point_type& p, const int level = -1) {
362:       this->_indexArray->clear();
363:       int size = 0;

365:       if (level == 0) {
366:         const index_type& idx = section->getIndex(p);

368:         this->_indexArray->push_back(oIndex_type(idx, 0));
369:         size += std::abs(idx.prefix);
370:       } else if ((level == 1) || (this->height() == 1)) {
371:         const Obj<typename sieve_type::coneSequence>& cone = this->_sieve->cone(p);
372:         typename sieve_type::coneSequence::iterator   end  = cone->end();
373:         const index_type& idx = section->getIndex(p);

375:         this->_indexArray->push_back(idx);
376:         size += idx.prefix;
377:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
378:           const index_type& pIdx = section->getIndex(*p_iter);

380:           this->_indexArray->push_back(oIndex_type(pIdx, 0));
381:           size += std::abs(pIdx.prefix);
382:         }
383:       } else if (level == -1) {
384:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
385:         typename oConeArray::iterator end     = closure->end();

387:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
388:           const index_type& pIdx = section->getIndex(p_iter->first);

390:           this->_indexArray->push_back(oIndex_type(pIdx, p_iter->second));
391:           size += std::abs(pIdx.prefix);
392:         }
393:       } else {
394:         throw ALE::Exception("Bundle has not yet implemented nCone");
395:       }
396:       if (this->debug()) {
397:         for(typename oIndexArray::iterator i_iter = this->_indexArray->begin(); i_iter != this->_indexArray->end(); ++i_iter) {
398:           printf("[%d]index interval (%d, %d)\n", this->commRank(), i_iter->first.prefix, i_iter->first.index);
399:         }
400:       }
401:       int *indexArray = this->getIndexArray(size);

403:       this->expandIntervals(this->_indexArray, indexArray);
404:       return indices_type(indexArray, size);
405:     };
406:     template<typename Section_, typename Numbering>
407:     const indices_type getIndices(const Obj<Section_>& section, const point_type& p, const Obj<Numbering>& numbering, const int level = -1) {
408:       int *indexArray = NULL;
409:       int  size       = 0;

411:       if (level == 0) {
412:         size      += section->getFiberDimension(p);
413:         indexArray = this->getIndexArray(size);
414:         int  k     = 0;

416:         section->getIndices(p, numbering, indexArray, &k);
417:       } else if ((level == 1) || (this->height() == 1)) {
418:         const Obj<typename sieve_type::coneSequence>& cone = this->_sieve->cone(p);
419:         typename sieve_type::coneSequence::iterator   end  = cone->end();

421:         size      += section->getFiberDimension(p);
422:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
423:           size    += section->getFiberDimension(*p_iter);
424:         }
425:         indexArray = this->getIndexArray(size);
426:         int  k     = 0;

428:         section->getIndices(p, numbering, indexArray, &k);
429:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
430:           section->getIndices(*p_iter, numbering, indexArray, &k);
431:         }
432:       } else if (level == -1) {
433:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
434:         typename oConeArray::iterator end     = closure->end();

436:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
437:           size    += section->getFiberDimension(p_iter->first);
438:         }
439:         indexArray = this->getIndexArray(size);
440:         int  k     = 0;
441:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
442:           section->getIndices(p_iter->first, numbering, indexArray, &k, p_iter->second);
443:         }
444:       } else {
445:         throw ALE::Exception("Bundle has not yet implemented getIndices() for an arbitrary level");
446:       }
447:       return indices_type(indexArray, size);
448:     };
449:   public: // Retrieval traversal
450:     // Return the values for the closure of this point
451:     //   use a smart pointer?
452:     template<typename Section_>
453:     const typename Section_::value_type *restrict(const Obj<Section_>& section, const point_type& p) {
454:       const int                       size   = this->sizeWithBC(section, p);
455:       typename Section_::value_type  *values = section->getRawArray(size);
456:       int                             j      = -1;

458:       // We could actually ask for the height of the individual point
459:       if (this->height() < 2) {
460:         const int& dim = section->getFiberDimension(p);
461:         const typename Section_::value_type *array = section->restrictPoint(p);

463:         for(int i = 0; i < dim; ++i) {
464:           values[++j] = array[i];
465:         }
466:         const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
467:         typename sieve_type::coneSequence::iterator   end  = cone->end();

469:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
470:           const int& dim = section->getFiberDimension(*p_iter);

472:           array = section->restrictPoint(*p_iter);
473:           for(int i = 0; i < dim; ++i) {
474:             values[++j] = array[i];
475:           }
476:         }
477:       } else {
478:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
479:         typename oConeArray::iterator end     = closure->end();

481:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
482:           const typename Section_::value_type *array = section->restrictPoint(p_iter->first);
483:           const int& dim = section->getFiberDimension(p_iter->first);

485:           if (p_iter->second >= 0) {
486:             for(int i = 0; i < dim; ++i) {
487:               values[++j] = array[i];
488:             }
489:           } else {
490:             for(int i = dim-1; i >= 0; --i) {
491:               values[++j] = array[i];
492:             }
493:           }
494:         }
495:       }
496:       if (j != size-1) {
497:         ostringstream txt;

499:         txt << "Invalid restrict to point " << p << std::endl;
500:         txt << "  j " << j << " should be " << (size-1) << std::endl;
501:         std::cout << txt.str();
502:         throw ALE::Exception(txt.str().c_str());
503:       }
504:       return values;
505:     };
506:     template<typename Section_>
507:     void update(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
508:       int j = 0;

510:       if (this->height() < 2) {
511:         section->updatePoint(p, &v[j]);
512:         j += section->getFiberDimension(p);
513:         const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);

515:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
516:           section->updatePoint(*p_iter, &v[j]);
517:           j += section->getFiberDimension(*p_iter);
518:         }
519:       } else {
520:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
521:         typename oConeArray::iterator end     = closure->end();

523:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
524:           section->updatePoint(p_iter->first, &v[j], p_iter->second);
525:           j += section->getFiberDimension(p_iter->first);
526:         }
527:       }
528:     };
529:     template<typename Section_>
530:     void updateAdd(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
531:       int j = 0;

533:       if (this->height() < 2) {
534:         section->updateAddPoint(p, &v[j]);
535:         j += section->getFiberDimension(p);
536:         const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);

538:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
539:           section->updateAddPoint(*p_iter, &v[j]);
540:           j += section->getFiberDimension(*p_iter);
541:         }
542:       } else {
543:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
544:         typename oConeArray::iterator end     = closure->end();

546:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
547:           section->updateAddPoint(p_iter->first, &v[j], p_iter->second);
548:           j += section->getFiberDimension(p_iter->first);
549:         }
550:       }
551:     };
552:     template<typename Section_>
553:     void updateBC(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
554:       int j = 0;

556:       if (this->height() < 2) {
557:         section->updatePointBC(p, &v[j]);
558:         j += section->getFiberDimension(p);
559:         const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);

561:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
562:           section->updatePointBC(*p_iter, &v[j]);
563:           j += section->getFiberDimension(*p_iter);
564:         }
565:       } else {
566:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
567:         typename oConeArray::iterator end     = closure->end();

569:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
570:           section->updatePointBC(p_iter->first, &v[j], p_iter->second);
571:           j += section->getFiberDimension(p_iter->first);
572:         }
573:       }
574:     };
575:     template<typename Section_>
576:     void updateAll(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
577:       int j = 0;

579:       if (this->height() < 2) {
580:         section->updatePointBC(p, &v[j]);
581:         j += section->getFiberDimension(p);
582:         const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);

584:         for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
585:           section->updatePointAll(*p_iter, &v[j]);
586:           j += section->getFiberDimension(*p_iter);
587:         }
588:       } else {
589:         const Obj<oConeArray>         closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
590:         typename oConeArray::iterator end     = closure->end();

592:         for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
593:           section->updatePointAll(p_iter->first, &v[j], p_iter->second);
594:           j += section->getFiberDimension(p_iter->first);
595:         }
596:       }
597:     };
598:   public: // Allocation
599:     template<typename Section_>
600:     void allocate(const Obj<Section_>& section, const Obj<send_overlap_type>& sendOverlap = NULL) {
601:       bool doGhosts = !sendOverlap.isNull();

603:       this->_factory->orderPatch(section, this->getSieve(), sendOverlap);
604:       if (doGhosts) {
605:         if (this->_debug > 1) {std::cout << "Ordering patch for ghosts" << std::endl;}
606:         const typename Section_::chart_type& points = section->getChart();
607:         int offset = 0;

609:         for(typename Section_::chart_type::iterator point = points.begin(); point != points.end(); ++point) {
610:           const typename Section_::index_type& idx = section->getIndex(*point);

612:           offset = std::max(offset, idx.index + std::abs(idx.prefix));
613:         }
614:         this->_factory->orderPatch(section, this->getSieve(), NULL, offset);
615:         if (offset != section->sizeWithBC()) throw ALE::Exception("Inconsistent array sizes in section");
616:       }
617:       section->allocateStorage();
618:     };
619:     template<typename Section_>
620:     void reallocate(const Obj<Section_>& section) {
621:       if (section->getNewAtlas().isNull()) return;
622:       // Since copy() preserves offsets, we must reinitialize them before ordering
623:       const Obj<typename Section_::atlas_type>         atlas    = section->getAtlas();
624:       const Obj<typename Section_::atlas_type>&        newAtlas = section->getNewAtlas();
625:       const typename Section_::atlas_type::chart_type& chart    = newAtlas->getChart();
626:       const typename Section_::atlas_type::chart_type& oldChart = atlas->getChart();
627:       int                                              newSize  = 0;
628:       typename Section_::index_type                    defaultIdx(0, -1);

630:       for(typename Section_::atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
631:         defaultIdx.prefix = newAtlas->restrictPoint(*c_iter)[0].prefix;
632:         newAtlas->updatePoint(*c_iter, &defaultIdx);
633:         newSize += defaultIdx.prefix;
634:       }
635:       section->setAtlas(newAtlas);
636:       this->_factory->orderPatch(section, this->getSieve());
637:       // Copy over existing values
638:       typename Section_::value_type                   *newArray = new typename Section_::value_type[newSize];
639:       const typename Section_::value_type             *array    = section->restrict();

641:       for(typename Section_::atlas_type::chart_type::const_iterator c_iter = oldChart.begin(); c_iter != oldChart.end(); ++c_iter) {
642:         const int& dim       = section->getFiberDimension(*c_iter);
643:         const int& offset    = atlas->restrictPoint(*c_iter)->index;
644:         const int& newOffset = newAtlas->restrictPoint(*c_iter)->index;

646:         for(int i = 0; i < dim; ++i) {
647:           newArray[newOffset+i] = array[offset+i];
648:         }
649:       }
650:       section->replaceStorage(newArray);
651:     };
652:   public: // Overlap
653:     template<typename Sequence>
654:     void constructOverlap(const Obj<Sequence>& points, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
655:       point_type *sendBuf = new point_type[points->size()];
656:       int         size    = 0;
657:       for(typename Sequence::iterator l_iter = points->begin(); l_iter != points->end(); ++l_iter) {
658:         sendBuf[size++] = *l_iter;
659:       }
660:       int *sizes   = new int[this->commSize()];   // The number of points coming from each process
661:       int *offsets = new int[this->commSize()+1]; // Prefix sums for sizes
662:       int *oldOffs = new int[this->commSize()+1]; // Temporary storage
663:       point_type *remotePoints = NULL;            // The points from each process
664:       int        *remoteRanks  = NULL;            // The rank and number of overlap points of each process that overlaps another

666:       // Change to Allgather() for the correct binning algorithm
667:       MPI_Gather(&size, 1, MPI_INT, sizes, 1, MPI_INT, 0, this->comm());
668:       if (this->commRank() == 0) {
669:         offsets[0] = 0;
670:         for(int p = 1; p <= this->commSize(); p++) {
671:           offsets[p] = offsets[p-1] + sizes[p-1];
672:         }
673:         remotePoints = new point_type[offsets[this->commSize()]];
674:       }
675:       MPI_Gatherv(sendBuf, size, MPI_INT, remotePoints, sizes, offsets, MPI_INT, 0, this->comm());
676:       std::map<int, std::map<int, std::set<point_type> > > overlapInfo; // Maps (p,q) to their set of overlap points

678:       if (this->commRank() == 0) {
679:         for(int p = 0; p < this->commSize(); p++) {
680:           std::sort(&remotePoints[offsets[p]], &remotePoints[offsets[p+1]]);
681:         }
682:         for(int p = 0; p <= this->commSize(); p++) {
683:           oldOffs[p] = offsets[p];
684:         }
685:         for(int p = 0; p < this->commSize(); p++) {
686:           for(int q = p+1; q < this->commSize(); q++) {
687:             std::set_intersection(&remotePoints[oldOffs[p]], &remotePoints[oldOffs[p+1]],
688:                                   &remotePoints[oldOffs[q]], &remotePoints[oldOffs[q+1]],
689:                                   std::insert_iterator<std::set<point_type> >(overlapInfo[p][q], overlapInfo[p][q].begin()));
690:             overlapInfo[q][p] = overlapInfo[p][q];
691:           }
692:           sizes[p]     = overlapInfo[p].size()*2;
693:           offsets[p+1] = offsets[p] + sizes[p];
694:         }
695:         remoteRanks = new int[offsets[this->commSize()]];
696:         int       k = 0;
697:         for(int p = 0; p < this->commSize(); p++) {
698:           for(typename std::map<int, std::set<point_type> >::iterator r_iter = overlapInfo[p].begin(); r_iter != overlapInfo[p].end(); ++r_iter) {
699:             remoteRanks[k*2]   = r_iter->first;
700:             remoteRanks[k*2+1] = r_iter->second.size();
701:             k++;
702:           }
703:         }
704:       }
705:       int numOverlaps;                          // The number of processes overlapping this process
706:       MPI_Scatter(sizes, 1, MPI_INT, &numOverlaps, 1, MPI_INT, 0, this->comm());
707:       int *overlapRanks = new int[numOverlaps]; // The rank and overlap size for each overlapping process
708:       MPI_Scatterv(remoteRanks, sizes, offsets, MPI_INT, overlapRanks, numOverlaps, MPI_INT, 0, this->comm());
709:       point_type *sendPoints = NULL;            // The points to send to each process
710:       if (this->commRank() == 0) {
711:         for(int p = 0, k = 0; p < this->commSize(); p++) {
712:           sizes[p] = 0;
713:           for(int r = 0; r < (int) overlapInfo[p].size(); r++) {
714:             sizes[p] += remoteRanks[k*2+1];
715:             k++;
716:           }
717:           offsets[p+1] = offsets[p] + sizes[p];
718:         }
719:         sendPoints = new point_type[offsets[this->commSize()]];
720:         for(int p = 0, k = 0; p < this->commSize(); p++) {
721:           for(typename std::map<int, std::set<point_type> >::iterator r_iter = overlapInfo[p].begin(); r_iter != overlapInfo[p].end(); ++r_iter) {
722:             int rank = r_iter->first;
723:             for(typename std::set<point_type>::iterator p_iter = (overlapInfo[p][rank]).begin(); p_iter != (overlapInfo[p][rank]).end(); ++p_iter) {
724:               sendPoints[k++] = *p_iter;
725:             }
726:           }
727:         }
728:       }
729:       int numOverlapPoints = 0;
730:       for(int r = 0; r < numOverlaps/2; r++) {
731:         numOverlapPoints += overlapRanks[r*2+1];
732:       }
733:       point_type *overlapPoints = new point_type[numOverlapPoints];
734:       MPI_Scatterv(sendPoints, sizes, offsets, MPI_INT, overlapPoints, numOverlapPoints, MPI_INT, 0, this->comm());

736:       for(int r = 0, k = 0; r < numOverlaps/2; r++) {
737:         int rank = overlapRanks[r*2];

739:         for(int p = 0; p < overlapRanks[r*2+1]; p++) {
740:           point_type point = overlapPoints[k++];

742:           sendOverlap->addArrow(point, rank, point);
743:           recvOverlap->addArrow(rank, point, point);
744:         }
745:       }

747:       delete [] overlapPoints;
748:       delete [] overlapRanks;
749:       delete [] sizes;
750:       delete [] offsets;
751:       delete [] oldOffs;
752:       if (this->commRank() == 0) {
753:         delete [] remoteRanks;
754:         delete [] remotePoints;
755:         delete [] sendPoints;
756:       }
757:     };
758:     void constructOverlap() {
759:       if (this->_calculatedOverlap) return;
760:       this->constructOverlap(this->getSieve()->base(), this->getSendOverlap(), this->getRecvOverlap());
761:       this->constructOverlap(this->getSieve()->cap(),  this->getSendOverlap(), this->getRecvOverlap());
762:       if (this->debug()) {
763:         this->_sendOverlap->view("Send overlap");
764:         this->_recvOverlap->view("Receive overlap");
765:       }
766:       this->_calculatedOverlap = true;
767:     };
768:   };
769:   class Discretization : public ALE::ParallelObject {
770:   public:
771:     typedef std::vector<Obj<Discretization> > children_type;
772:   protected:
773:     children_type _children;
774:     std::map<int,int> _dim2dof;
775:     std::map<int,int> _dim2class;
776:     const double *_points;
777:     const double *_weights;
778:     const double *_basis;
779:     const double *_basisDer;
780:     const int    *_indices;
781:   public:
782:     Discretization(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _points(NULL), _weights(NULL), _basis(NULL), _basisDer(NULL), _indices(NULL) {};
783:     virtual ~Discretization() {
784:       if (this->_indices) {delete [] this->_indices;}
785:     };
786:   public:
787:     const double *getQuadraturePoints() {return this->_points;};
788:     void          setQuadraturePoints(const double *points) {this->_points = points;};
789:     const double *getQuadratureWeights() {return this->_weights;};
790:     void          setQuadratureWeights(const double *weights) {this->_weights = weights;};
791:     const double *getBasis() {return this->_basis;};
792:     void          setBasis(const double *basis) {this->_basis = basis;};
793:     const double *getBasisDerivatives() {return this->_basisDer;};
794:     void          setBasisDerivatives(const double *basisDer) {this->_basisDer = basisDer;};
795:     int  getNumDof(const int dim) {return this->_dim2dof[dim];};
796:     void setNumDof(const int dim, const int numDof) {this->_dim2dof[dim] = numDof;};
797:     int  getDofClass(const int dim) {return this->_dim2class[dim];};
798:     void setDofClass(const int dim, const int dofClass) {this->_dim2class[dim] = dofClass;};
799:   public:
800:     void addChild(const Obj<Discretization>& child) {this->_children.push_back(child);};
801:     children_type getDiscs() {
802:       children_type       discs;
803:       Obj<Discretization> me = this;
804:       me.addRef();

806:       discs.push_back(me);
807:       for(children_type::iterator c_iter = this->_children.begin(); c_iter != this->_children.end(); ++c_iter) {
808:         const children_type childDiscs = (*c_iter)->getDiscs();
809:         discs.insert(discs.end(), childDiscs.begin(), childDiscs.end());
810:       }
811:       return discs;
812:     };
813:     const int *getIndices() {return this->_indices;};
814:     void       setIndices(const int *indices) {this->_indices = indices;};
815:     template<typename Bundle>
816:     int size(const Obj<Bundle>& mesh) {
817:       const Obj<typename Bundle::label_sequence>& cells   = mesh->heightStratum(0);
818:       const Obj<typename Bundle::coneArray>       closure = ALE::SieveAlg<Bundle>::closure(mesh, *cells->begin());
819:       const typename Bundle::coneArray::iterator  end     = closure->end();
820:       int                                         size    = 0;

822:       for(typename Bundle::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
823:         size += this->_dim2dof[mesh->depth(*cl_iter)];
824:       }
825:       return size;
826:     };
827:     template<typename Bundle>
828:     void calculateIndices(const Obj<Bundle>& mesh) {
829:       // Should have an iterator over the whole tree
830:       children_type discs = this->getDiscs();
831:       std::map<Discretization*, std::pair<int, int*> > indices;

833:       for(typename children_type::iterator c_iter = discs.begin(); c_iter != discs.end(); ++c_iter) {
834:         indices[*c_iter] = std::pair<int, int*>(0, new int[(*c_iter)->size(mesh)]);
835:         (*c_iter)->setIndices(indices[*c_iter].second);
836:       }
837:       const Obj<typename Bundle::label_sequence>& cells   = mesh->heightStratum(0);
838:       const Obj<typename Bundle::coneArray>       closure = ALE::SieveAlg<Bundle>::closure(mesh, *cells->begin());
839:       const typename Bundle::coneArray::iterator  end     = closure->end();
840:       int                                         offset  = 0;

842:       std::cout << "Closure for first element" << std::endl;
843:       for(typename Bundle::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
844:         const int dim = mesh->depth(*cl_iter);

846:         std::cout << "  point " << *cl_iter << " depth " << dim << std::endl;
847:         for(typename children_type::iterator c_iter = discs.begin(); c_iter != discs.end(); ++c_iter) {
848:           const int num = (*c_iter)->getNumDof(dim);

850:           std::cout << "    disc " << (*c_iter)->getName() << " numDof " << num << std::endl;
851:           for(int o = 0; o < num; ++o) {
852:             indices[*c_iter].second[indices[*c_iter].first++] = offset++;
853:           }
854:         }
855:       }
856:       for(typename children_type::iterator c_iter = discs.begin(); c_iter != discs.end(); ++c_iter) {
857:         std::cout << "Discretization " << (*c_iter)->getName() << " indices:";
858:         for(int i = 0; i < indices[*c_iter].first; ++i) {
859:           std::cout << " " << indices[*c_iter].second[i];
860:         }
861:         std::cout << std::endl;
862:       }
863:     };
864:   };
865:   class BoundaryCondition : public ALE::ParallelObject {
866:   public:
867:     typedef double (*function_type)(const double []);
868:     typedef double (*integrator_type)(const double [], const double [], const int, function_type);
869:   protected:
870:     std::string     _labelName;
871:     function_type   _func;
872:     integrator_type _integrator;
873:   public:
874:     BoundaryCondition(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
875:     ~BoundaryCondition() {};
876:   public:
877:     const std::string& getLabelName() {return this->_labelName;};
878:     void setLabelName(const std::string& name) {this->_labelName = name;};
879:     void setFunction(function_type func) {this->_func = func;};
880:     integrator_type getDualIntegrator() {return this->_integrator;};
881:     void setDualIntegrator(integrator_type integrator) {this->_integrator = integrator;};
882:   public:
883:     double evaluate(const double coords[]) {return this->_func(coords);};
884:     double integrateDual(const double v0[], const double J[], const int dualIndex) {return this->_integrator(v0, J, dualIndex, this->_func);};
885:   };
886:   class DiscretizationNew : public ALE::ParallelObject {
887:   protected:
888:     Obj<BoundaryCondition> _boundaryCondition;
889:     Obj<BoundaryCondition> _exactSolution;
890:     std::map<int,int> _dim2dof;
891:     std::map<int,int> _dim2class;
892:     const double *_points;
893:     const double *_weights;
894:     const double *_basis;
895:     const double *_basisDer;
896:     const int    *_indices;
897:   public:
898:     DiscretizationNew(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _points(NULL), _weights(NULL), _basis(NULL), _basisDer(NULL), _indices(NULL) {};
899:     virtual ~DiscretizationNew() {
900:       if (this->_indices) {delete [] this->_indices;}
901:     };
902:   public:
903:     const Obj<BoundaryCondition>& getBoundaryCondition() {return this->_boundaryCondition;};
904:     void setBoundaryCondition(const Obj<BoundaryCondition>& boundaryCondition) {this->_boundaryCondition = boundaryCondition;};
905:     const Obj<BoundaryCondition>& getExactSolution() {return this->_exactSolution;};
906:     void setExactSolution(const Obj<BoundaryCondition>& exactSolution) {this->_exactSolution = exactSolution;};
907:     const double *getQuadraturePoints() {return this->_points;};
908:     void          setQuadraturePoints(const double *points) {this->_points = points;};
909:     const double *getQuadratureWeights() {return this->_weights;};
910:     void          setQuadratureWeights(const double *weights) {this->_weights = weights;};
911:     const double *getBasis() {return this->_basis;};
912:     void          setBasis(const double *basis) {this->_basis = basis;};
913:     const double *getBasisDerivatives() {return this->_basisDer;};
914:     void          setBasisDerivatives(const double *basisDer) {this->_basisDer = basisDer;};
915:     int  getNumDof(const int dim) {return this->_dim2dof[dim];};
916:     void setNumDof(const int dim, const int numDof) {this->_dim2dof[dim] = numDof;};
917:     int  getDofClass(const int dim) {return this->_dim2class[dim];};
918:     void setDofClass(const int dim, const int dofClass) {this->_dim2class[dim] = dofClass;};
919:   public:
920:     const int *getIndices() {return this->_indices;};
921:     void       setIndices(const int *indices) {this->_indices = indices;};
922:     template<typename Bundle>
923:     int size(const Obj<Bundle>& mesh) {
924:       const Obj<typename Bundle::label_sequence>& cells   = mesh->heightStratum(0);
925:       const Obj<typename Bundle::coneArray>       closure = ALE::SieveAlg<Bundle>::closure(mesh, *cells->begin());
926:       const typename Bundle::coneArray::iterator  end     = closure->end();
927:       int                                         size    = 0;

929:       for(typename Bundle::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
930:         size += this->_dim2dof[mesh->depth(*cl_iter)];
931:       }
932:       return size;
933:     };
934:   };
935: }

937: namespace ALE {
938: #ifdef NEW_SECTION
939:   class Mesh : public Bundle<ALE::Sieve<int,int,int>, GeneralSection<int, double> > {
940: #else
941:   class Mesh : public Bundle<ALE::Sieve<int,int,int> > {
942: #endif
943:   public:
944: #ifdef NEW_SECTION
945:     typedef Bundle<ALE::Sieve<int,int,int>, GeneralSection<int, double> > base_type;
946: #else
947:     typedef Bundle<ALE::Sieve<int,int,int> > base_type;
948: #endif
949:     typedef base_type::sieve_type            sieve_type;
950:     typedef sieve_type::point_type           point_type;
951:     typedef base_type::label_sequence        label_sequence;
952:     typedef base_type::real_section_type     real_section_type;
953:     typedef base_type::numbering_type        numbering_type;
954:     typedef base_type::order_type            order_type;
955:     typedef base_type::send_overlap_type     send_overlap_type;
956:     typedef base_type::recv_overlap_type     recv_overlap_type;
957:     typedef base_type::sieve_alg_type        sieve_alg_type;
958:     typedef std::map<std::string, Obj<DiscretizationNew> > discretizations_type;
959:   protected:
960:     int                   _dim;
961:     // Discretization
962:     Obj<Discretization>    _discretization;
963:     Obj<BoundaryCondition> _boundaryCondition;
964:     discretizations_type   _discretizations;
965:   public:
966:     Mesh(MPI_Comm comm, int dim, int debug = 0) : base_type(comm, debug), _dim(dim) {
967:       ///this->_factory = NumberingFactory::singleton(debug);
968:       this->_discretization = new Discretization(comm, debug);
969:       this->_boundaryCondition = new BoundaryCondition(comm, debug);
970:     };
971:   public: // Accessors
972:     int getDimension() const {return this->_dim;};
973:     void setDimension(const int dim) {this->_dim = dim;};
974:     const Obj<Discretization>&    getDiscretization() {return this->_discretization;};
975:     void setDiscretization(const Obj<Discretization>& discretization) {this->_discretization = discretization;};
976:     const Obj<BoundaryCondition>& getBoundaryCondition() {return this->_boundaryCondition;};
977:     void setBoundaryCondition(const Obj<BoundaryCondition>& boundaryCondition) {this->_boundaryCondition = boundaryCondition;};
978:     const Obj<DiscretizationNew>& getDiscretization(const std::string& name) {return this->_discretizations[name];};
979:     void setDiscretization(const std::string& name, const Obj<DiscretizationNew>& disc) {this->_discretizations[name] = disc;};
980:     Obj<std::set<std::string> > getDiscretizations() const {
981:       Obj<std::set<std::string> > names = std::set<std::string>();

983:       for(discretizations_type::const_iterator d_iter = this->_discretizations.begin(); d_iter != this->_discretizations.end(); ++d_iter) {
984:         names->insert(d_iter->first);
985:       }
986:       return names;
987:     };
988:   public: // Mesh geometry
989:     void computeTriangleGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
990:       const double    *coords = this->restrict(coordinates, e);
991:       const int        dim    = 2;
992:       double           invDet;

994:       if (v0) {
995:         for(int d = 0; d < dim; d++) {
996:           v0[d] = coords[d];
997:         }
998:       }
999:       if (J) {
1000:         for(int d = 0; d < dim; d++) {
1001:           for(int f = 0; f < dim; f++) {
1002:             J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
1003:           }
1004:         }
1005:         detJ = J[0]*J[3] - J[1]*J[2];
1006:       }
1007:       if (invJ) {
1008:         invDet  = 1.0/detJ;
1009:         invJ[0] =  invDet*J[3];
1010:         invJ[1] = -invDet*J[1];
1011:         invJ[2] = -invDet*J[2];
1012:         invJ[3] =  invDet*J[0];
1013:       }
1014:     };
1015:     void computeQuadrilateralGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double point[], double v0[], double J[], double invJ[], double& detJ) {
1016:       const double *coords = this->restrict(coordinates, e);
1017:       const int     dim    = 2;
1018:       double        invDet;

1020:       if (v0) {
1021:         for(int d = 0; d < dim; d++) {
1022:           v0[d] = coords[d];
1023:         }
1024:       }
1025:       if (J) {
1026:         double x_1 = coords[2] - coords[0];
1027:         double y_1 = coords[3] - coords[1];
1028:         double x_2 = coords[6] - coords[0];
1029:         double y_2 = coords[7] - coords[1];
1030:         double x_3 = coords[4] - coords[0];
1031:         double y_3 = coords[5] - coords[1];

1033:         J[0] = x_1 + (x_3 - x_1 - x_2)*point[1];
1034:         J[1] = x_2 + (x_3 - x_1 - x_2)*point[0];
1035:         J[2] = y_1 + (y_3 - y_1 - y_2)*point[1];
1036:         J[3] = y_1 + (y_3 - y_1 - y_2)*point[0];
1037:         detJ = J[0]*J[3] - J[1]*J[2];
1038:       }
1039:       if (invJ) {
1040:         invDet  = 1.0/detJ;
1041:         invJ[0] =  invDet*J[3];
1042:         invJ[1] = -invDet*J[1];
1043:         invJ[2] = -invDet*J[2];
1044:         invJ[3] =  invDet*J[0];
1045:       }
1046:     };
1047:     void computeTetrahedronGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1048:       const double *coords = this->restrict(coordinates, e);
1049:       const int     dim    = 3;
1050:       double        invDet;

1052:       if (v0) {
1053:         for(int d = 0; d < dim; d++) {
1054:           v0[d] = coords[d];
1055:         }
1056:       }
1057:       if (J) {
1058:         for(int d = 0; d < dim; d++) {
1059:           for(int f = 0; f < dim; f++) {
1060:             J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
1061:           }
1062:         }
1063:         // The minus sign is here since I orient the first face to get the outward normal
1064:         detJ = -(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
1065:                  J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
1066:                  J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
1067:       }
1068:       if (invJ) {
1069:         invDet  = -1.0/detJ;
1070:          invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
1071:          invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
1072:          invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
1073:          invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
1074:          invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
1075:          invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
1076:          invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
1077:          invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
1078:          invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
1079:       }
1080:     };
1081:     void computeHexahedralGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double point[], double v0[], double J[], double invJ[], double& detJ) {
1082:       const double *coords = this->restrict(coordinates, e);
1083:       const int     dim    = 3;
1084:       double        invDet;

1086:       if (v0) {
1087:         for(int d = 0; d < dim; d++) {
1088:           v0[d] = coords[d];
1089:         }
1090:       }
1091:       if (J) {
1092:         double x_1 = coords[3]  - coords[0];
1093:         double y_1 = coords[4]  - coords[1];
1094:         double z_1 = coords[5]  - coords[2];
1095:         double x_2 = coords[6]  - coords[0];
1096:         double y_2 = coords[7]  - coords[1];
1097:         double z_2 = coords[8]  - coords[2];
1098:         double x_3 = coords[9]  - coords[0];
1099:         double y_3 = coords[10] - coords[1];
1100:         double z_3 = coords[11] - coords[2];
1101:         double x_4 = coords[12] - coords[0];
1102:         double y_4 = coords[13] - coords[1];
1103:         double z_4 = coords[14] - coords[2];
1104:         double x_5 = coords[15] - coords[0];
1105:         double y_5 = coords[16] - coords[1];
1106:         double z_5 = coords[17] - coords[2];
1107:         double x_6 = coords[18] - coords[0];
1108:         double y_6 = coords[19] - coords[1];
1109:         double z_6 = coords[20] - coords[2];
1110:         double x_7 = coords[21] - coords[0];
1111:         double y_7 = coords[22] - coords[1];
1112:         double z_7 = coords[23] - coords[2];
1113:         double g_x = x_1 - x_2 + x_3 + x_4 - x_5 + x_6 - x_7;
1114:         double g_y = y_1 - y_2 + y_3 + y_4 - y_5 + y_6 - y_7;
1115:         double g_z = z_1 - z_2 + z_3 + z_4 - z_5 + z_6 - z_7;

1117:         J[0] = x_1 + (x_2 - x_1 - x_3)*point[1] + (x_5 - x_1 - x_4)*point[2] + g_x*point[1]*point[2];
1118:         J[1] = x_3 + (x_2 - x_1 - x_3)*point[0] + (x_7 - x_3 - x_4)*point[2] + g_x*point[2]*point[0];
1119:         J[2] = x_4 + (x_7 - x_3 - x_4)*point[1] + (x_5 - x_1 - x_4)*point[0] + g_x*point[0]*point[1];
1120:         J[3] = y_1 + (y_2 - y_1 - y_3)*point[1] + (y_5 - y_1 - y_4)*point[2] + g_y*point[1]*point[2];
1121:         J[4] = y_3 + (y_2 - y_1 - y_3)*point[0] + (y_7 - y_3 - y_4)*point[2] + g_y*point[2]*point[0];
1122:         J[5] = y_4 + (y_7 - y_3 - y_4)*point[1] + (y_5 - y_1 - y_4)*point[0] + g_y*point[0]*point[1];
1123:         J[6] = z_1 + (z_2 - z_1 - z_3)*point[1] + (z_5 - z_1 - z_4)*point[2] + g_z*point[1]*point[2];
1124:         J[7] = z_3 + (z_2 - z_1 - z_3)*point[0] + (z_7 - z_3 - z_4)*point[2] + g_z*point[2]*point[0];
1125:         J[8] = z_4 + (z_7 - z_3 - z_4)*point[1] + (z_5 - z_1 - z_4)*point[0] + g_z*point[0]*point[1];
1126:         detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
1127:                 J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
1128:                 J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
1129:       }
1130:       if (invJ) {
1131:         invDet  = 1.0/detJ;
1132:          invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
1133:          invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
1134:          invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
1135:          invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
1136:          invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
1137:          invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
1138:          invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
1139:          invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
1140:          invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
1141:       }
1142:     };
1143:     void computeElementGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1144:       if (this->_dim == 2) {
1145:         computeTriangleGeometry(coordinates, e, v0, J, invJ, detJ);
1146:       } else if (this->_dim == 3) {
1147:         computeTetrahedronGeometry(coordinates, e, v0, J, invJ, detJ);
1148:       } else {
1149:         throw ALE::Exception("Unsupport dimension for element geometry computation");
1150:       }
1151:     };
1152:     double getMaxVolume() {
1153:       const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1154:       const Obj<label_sequence>&    cells       = this->heightStratum(0);
1155:       const int                     dim         = this->getDimension();
1156:       double v0[3], J[9], invJ[9], detJ, refVolume = 0.0, maxVolume = 0.0;

1158:       if (dim == 1) refVolume = 2.0;
1159:       if (dim == 2) refVolume = 2.0;
1160:       if (dim == 3) refVolume = 4.0/3.0;
1161:       for(label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
1162:         this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
1163:         maxVolume = std::max(maxVolume, detJ*refVolume);
1164:       }
1165:       return maxVolume;
1166:     };
1167:     // Find the cell in which this point lies (stupid algorithm)
1168:     point_type locatePoint_2D(const real_section_type::value_type point[]) {
1169:       const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1170:       const Obj<label_sequence>&    cells       = this->heightStratum(0);
1171:       const int                     embedDim    = 2;
1172:       double v0[2], J[4], invJ[4], detJ;

1174:       for(label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
1175:         this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
1176:         double xi   = invJ[0*embedDim+0]*(point[0] - v0[0]) + invJ[0*embedDim+1]*(point[1] - v0[1]);
1177:         double eta  = invJ[1*embedDim+0]*(point[0] - v0[0]) + invJ[1*embedDim+1]*(point[1] - v0[1]);

1179:         if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) {
1180:           return *c_iter;
1181:         }
1182:       }
1183:       throw ALE::Exception("Could not locate point");
1184:     };
1185:     //   Assume a simplex and 3D
1186:     point_type locatePoint_3D(const real_section_type::value_type point[]) {
1187:       const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1188:       const Obj<label_sequence>&    cells       = this->heightStratum(0);
1189:       const int                     embedDim    = 3;
1190:       double v0[3], J[9], invJ[9], detJ;

1192:       for(label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
1193:         this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
1194:         double xi   = invJ[0*embedDim+0]*(point[0] - v0[0]) + invJ[0*embedDim+1]*(point[1] - v0[1]) + invJ[0*embedDim+2]*(point[2] - v0[2]);
1195:         double eta  = invJ[1*embedDim+0]*(point[0] - v0[0]) + invJ[1*embedDim+1]*(point[1] - v0[1]) + invJ[1*embedDim+2]*(point[2] - v0[2]);
1196:         double zeta = invJ[2*embedDim+0]*(point[0] - v0[0]) + invJ[2*embedDim+1]*(point[1] - v0[1]) + invJ[2*embedDim+2]*(point[2] - v0[2]);

1198:         if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) {
1199:           return *c_iter;
1200:         }
1201:       }
1202:       throw ALE::Exception("Could not locate point");
1203:     };
1204:     point_type locatePoint(const real_section_type::value_type point[], point_type guess = -1) {
1205:       //guess overrides this by saying that we already know the relation of this point to this mesh.  We will need to make it a more robust "guess" later for more than P1
1206:       if (guess != -1) {
1207:         return guess;
1208:       }else if (this->_dim == 2) {
1209:         return locatePoint_2D(point);
1210:       } else if (this->_dim == 3) {
1211:         return locatePoint_3D(point);
1212:       } else {
1213:         throw ALE::Exception("No point location for mesh dimension");
1214:       }
1215:     };
1216:   public: // Discretization
1217:     void markBoundaryCells(const std::string& name) {
1218:       const Obj<label_type>&     label    = this->getLabel(name);
1219:       const Obj<label_sequence>& boundary = this->getLabelStratum(name, 1);
1220:       const Obj<sieve_type>&     sieve    = this->getSieve();

1222:       for(label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
1223:         if (this->height(*e_iter) == 1) {
1224:           const point_type cell = *sieve->support(*e_iter)->begin();

1226:           this->setValue(label, cell, 2);
1227:         }
1228:       }
1229:     };
1230:     void calculateIndices() {
1231:       // Should have an iterator over the whole tree
1232:       Obj<std::set<std::string> > discs = this->getDiscretizations();
1233:       Obj<Mesh>                   mesh  = this;
1234:       std::map<std::string, std::pair<int, int*> > indices;

1236:       mesh.addRef();
1237:       for(std::set<std::string>::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1238:         const Obj<DiscretizationNew>& disc = this->getDiscretization(*d_iter);

1240:         indices[*d_iter] = std::pair<int, int*>(0, new int[disc->size(mesh)]);
1241:         disc->setIndices(indices[*d_iter].second);
1242:       }
1243:       const Obj<label_sequence>& cells   = this->heightStratum(0);
1244:       const Obj<coneArray>       closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), *cells->begin());
1245:       const coneArray::iterator  end     = closure->end();
1246:       int                                         offset  = 0;

1248:       std::cout << "Closure for first element" << std::endl;
1249:       for(coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
1250:         const int dim = this->depth(*cl_iter);

1252:         std::cout << "  point " << *cl_iter << " depth " << dim << std::endl;
1253:         for(std::set<std::string>::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1254:           const Obj<DiscretizationNew>& disc = this->getDiscretization(*d_iter);
1255:           const int                  num  = disc->getNumDof(dim);

1257:           std::cout << "    disc " << disc->getName() << " numDof " << num << std::endl;
1258:           for(int o = 0; o < num; ++o) {
1259:             indices[*d_iter].second[indices[*d_iter].first++] = offset++;
1260:           }
1261:         }
1262:       }
1263:       for(std::set<std::string>::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1264:         const Obj<DiscretizationNew>& disc = this->getDiscretization(*d_iter);

1266:         std::cout << "Discretization " << disc->getName() << " indices:";
1267:         for(int i = 0; i < indices[*d_iter].first; ++i) {
1268:           std::cout << " " << indices[*d_iter].second[i];
1269:         }
1270:         std::cout << std::endl;
1271:       }
1272:     };
1273:     void setupField(const Obj<real_section_type>& s, const bool postponeGhosts = false) {
1274:       const std::string& name = this->_boundaryCondition->getLabelName();

1276:       for(int d = 0; d <= this->_dim; ++d) {
1277:         s->setFiberDimension(this->depthStratum(d), this->_discretization->getNumDof(d));
1278:       }
1279:       if (!name.empty()) {
1280:         const Obj<label_sequence>& boundary = this->getLabelStratum(name, 1);

1282:         for(label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
1283:           s->setConstraintDimension(*e_iter, this->_discretization->getNumDof(this->depth(*e_iter)));
1284:         }
1285:       }
1286:       if (postponeGhosts) throw ALE::Exception("Not implemented yet");
1287:       this->allocate(s);
1288:       s->defaultConstraintDof();
1289:       if (!name.empty()) {
1290:         const Obj<label_sequence>&     boundaryCells = this->getLabelStratum(name, 2);
1291:         const Obj<real_section_type>&  coordinates   = this->getRealSection("coordinates");
1292:         real_section_type::value_type *values        = new real_section_type::value_type[this->sizeWithBC(s, *boundaryCells->begin())];
1293:         double                        *v0            = new double[this->getDimension()];
1294:         double                        *J             = new double[this->getDimension()*this->getDimension()];
1295:         double                         detJ;

1297:         for(label_sequence::iterator c_iter = boundaryCells->begin(); c_iter != boundaryCells->end(); ++c_iter) {
1298:           const Obj<coneArray>      closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), *c_iter);
1299:           const coneArray::iterator end     = closure->end();
1300:           int                       v       = 0;

1302:           this->computeElementGeometry(coordinates, *c_iter, v0, J, PETSC_NULL, detJ);
1303:           for(coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
1304:             const int cDim = s->getConstraintDimension(*cl_iter);

1306:             if (cDim) {
1307:               for(int d = 0; d < cDim; ++d, ++v) {
1308:                 values[v] = this->_boundaryCondition->integrateDual(v0, J, v);
1309:               }
1310:             } else {
1311:               const int dim = s->getFiberDimension(*cl_iter);

1313:               for(int d = 0; d < dim; ++d, ++v) {
1314:                 values[v] = 0.0;
1315:               }
1316:             }
1317:           }
1318:           this->updateAll(s, *c_iter, values);
1319:         }
1320:         delete [] values;
1321:         delete [] v0;
1322:         delete [] J;
1323:       }
1324:     };
1325:     void setupFieldMultiple(const Obj<real_section_type>& s, const bool postponeGhosts = false) {
1326:       const Obj<std::set<std::string> >& discs = this->getDiscretizations();
1327:       std::set<std::string> names;

1329:       for(int d = 0; d <= this->_dim; ++d) {
1330:         int numDof = 0;

1332:         for(std::set<std::string>::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter) {
1333:           const Obj<ALE::DiscretizationNew>& disc = this->getDiscretization(*f_iter);
1334:           const Obj<ALE::BoundaryCondition>& bc   = disc->getBoundaryCondition();

1336:           numDof += disc->getNumDof(d);
1337:           if (!bc.isNull()) {
1338:             const Obj<label_sequence>& boundary = this->getLabelStratum(bc->getLabelName(), 1);

1340:             names.insert(bc->getLabelName());
1341:             for(label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
1342:               //s->addConstraintDimension(*e_iter, disc->getNumDof(this->depth(*e_iter)));
1343:             }
1344:           }
1345:         }
1346:         s->setFiberDimension(this->depthStratum(d), numDof);
1347:       }
1348:       if (postponeGhosts) throw ALE::Exception("Not implemented yet");
1349:       this->allocate(s);
1350:       s->defaultConstraintDof();
1351:       for(std::set<std::string>::const_iterator n_iter = names.begin(); n_iter != names.end(); ++n_iter) {
1352:         const Obj<label_sequence>&     boundaryCells = this->getLabelStratum(*n_iter, 2);
1353:         const Obj<real_section_type>&  coordinates   = this->getRealSection("coordinates");
1354:         real_section_type::value_type *values        = new real_section_type::value_type[this->sizeWithBC(s, *boundaryCells->begin())];
1355:         double                        *v0            = new double[this->getDimension()];
1356:         double                        *J             = new double[this->getDimension()*this->getDimension()];
1357:         double                         detJ;

1359:         for(label_sequence::iterator c_iter = boundaryCells->begin(); c_iter != boundaryCells->end(); ++c_iter) {
1360:           const Obj<coneArray>      closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), *c_iter);
1361:           const coneArray::iterator end     = closure->end();
1362:           int                       v       = 0;

1364:           this->computeElementGeometry(coordinates, *c_iter, v0, J, PETSC_NULL, detJ);
1365:           for(coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
1366:             const int cDim = s->getConstraintDimension(*cl_iter);

1368:             if (cDim) {
1369:               for(int d = 0; d < cDim; ++d, ++v) {
1370:                 values[v] = this->_boundaryCondition->integrateDual(v0, J, v);
1371:               }
1372:             } else {
1373:               const int dim = s->getFiberDimension(*cl_iter);

1375:               for(int d = 0; d < dim; ++d, ++v) {
1376:                 values[v] = 0.0;
1377:               }
1378:             }
1379:           }
1380:           this->updateAll(s, *c_iter, values);
1381:         }
1382:         delete [] values;
1383:         delete [] v0;
1384:         delete [] J;
1385:       }
1386:     };
1387:   public:
1388:     void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
1389:       if (comm == MPI_COMM_NULL) {
1390:         comm = this->comm();
1391:       }
1392:       if (name == "") {
1393:         PetscPrintf(comm, "viewing a Mesh\n");
1394:       } else {
1395:         PetscPrintf(comm, "viewing Mesh '%s'\n", name.c_str());
1396:       }
1397:       this->getSieve()->view("mesh sieve", comm);
1398:       Obj<std::set<std::string> > sections = this->getRealSections();

1400:       for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
1401:         this->getRealSection(*name)->view(*name);
1402:       }
1403:       sections = this->getIntSections();
1404:       for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
1405:         this->getIntSection(*name)->view(*name);
1406:       }
1407:       sections = this->getArrowSections();
1408:       for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
1409:         this->getArrowSection(*name)->view(*name);
1410:       }
1411:     };
1412:     template<typename value_type>
1413:     static std::string printMatrix(const std::string& name, const int rows, const int cols, const value_type matrix[], const int rank = -1)
1414:     {
1415:       ostringstream output;
1416:       ostringstream rankStr;

1418:       if (rank >= 0) {
1419:         rankStr << "[" << rank << "]";
1420:       }
1421:       output << rankStr.str() << name << " = " << std::endl;
1422:       for(int r = 0; r < rows; r++) {
1423:         if (r == 0) {
1424:           output << rankStr.str() << " /";
1425:         } else if (r == rows-1) {
1426:           output << rankStr.str() << " \\";
1427:         } else {
1428:           output << rankStr.str() << " |";
1429:         }
1430:         for(int c = 0; c < cols; c++) {
1431:           output << " " << matrix[r*cols+c];
1432:         }
1433:         if (r == 0) {
1434:           output << " \\" << std::endl;
1435:         } else if (r == rows-1) {
1436:           output << " /" << std::endl;
1437:         } else {
1438:           output << " |" << std::endl;
1439:         }
1440:       }
1441:       return output.str();
1442:     };
1443:   };
1444:   class MeshOld : public ALE::ParallelObject {
1445:   public:
1446:     // PCICE BC
1447:     typedef struct {double rho,u,v,p;}   bc_value_type;
1448:     typedef std::map<int, bc_value_type> bc_values_type;
1449:   protected:
1450:     int            _dim;
1451:     // PCICE BC
1452:     bc_values_type _bcValues;
1453:   public:
1454:     MeshOld(MPI_Comm comm, int dim, int debug = 0) : ALE::ParallelObject(comm, debug), _dim(dim) {};
1455:   public: // Accessors
1456: #if 0
1457:     // Only works in 2D
1458:     int orientation(const patch_type& patch, const point_type& cell) {
1459:       const Obj<real_section_type>&     coordinates = this->getRealSection("coordinates");
1460:       const Obj<topology_type>&         topology    = this->getTopology();
1461:       const Obj<sieve_type>&            sieve       = topology->getPatch(patch);
1462:       const Obj<sieve_type::coneArray>& cone        = sieve->nCone(cell, topology->depth());
1463:       sieve_type::coneArray::iterator   cBegin      = cone->begin();
1464:       real_section_type::value_type     root[2];
1465:       real_section_type::value_type     vA[2];
1466:       real_section_type::value_type     vB[2];

1468:       const real_section_type::value_type *coords = coordinates->restrictPoint(patch, *cBegin);
1469:       root[0] = coords[0];
1470:       root[1] = coords[1];
1471:       ++cBegin;
1472:       coords = coordinates->restrictPoint(patch, *cBegin);
1473:       vA[0] = coords[0] - root[0];
1474:       vA[1] = coords[1] - root[1];
1475:       ++cBegin;
1476:       coords = coordinates->restrictPoint(patch, *cBegin);
1477:       vB[0] = coords[0] - root[0];
1478:       vB[1] = coords[1] - root[1];
1479:       double det = vA[0]*vB[1] - vA[1]*vB[0];
1480:       if (det > 0.0) return  1;
1481:       if (det < 0.0) return -1;
1482:       return 0;
1483:     };
1484: #endif
1485:   public: // BC values for PCICE
1486:     const bc_value_type& getBCValue(const int bcFunc) {
1487:       return this->_bcValues[bcFunc];
1488:     };
1489:     void setBCValue(const int bcFunc, const bc_value_type& value) {
1490:       this->_bcValues[bcFunc] = value;
1491:     };
1492:     bc_values_type& getBCValues() {
1493:       return this->_bcValues;
1494:     };
1495:     void distributeBCValues() {
1496:       int size = this->_bcValues.size();

1498:       MPI_Bcast(&size, 1, MPI_INT, 0, this->comm());
1499:       if (this->commRank()) {
1500:         for(int bc = 0; bc < size; ++bc) {
1501:           int           funcNum;
1502:           bc_value_type funcVal;

1504:           MPI_Bcast((void *) &funcNum, 1, MPI_INT,    0, this->comm());
1505:           MPI_Bcast((void *) &funcVal, 4, MPI_DOUBLE, 0, this->comm());
1506:           this->_bcValues[funcNum] = funcVal;
1507:         }
1508:       } else {
1509:         for(bc_values_type::iterator bc_iter = this->_bcValues.begin(); bc_iter != this->_bcValues.end(); ++bc_iter) {
1510:           const int&           funcNum = bc_iter->first;
1511:           const bc_value_type& funcVal = bc_iter->second;
1512:           MPI_Bcast((void *) &funcNum, 1, MPI_INT,    0, this->comm());
1513:           MPI_Bcast((void *) &funcVal, 4, MPI_DOUBLE, 0, this->comm());
1514:         }
1515:       }
1516:     };
1517:   };
1518:   class MeshBuilder {
1519:     typedef ALE::Mesh Mesh;
1520:   public:
1523:     /*
1524:       Simple square boundary:

1526:      18--5-17--4--16
1527:       |     |     |
1528:       6    10     3
1529:       |     |     |
1530:      19-11-20--9--15
1531:       |     |     |
1532:       7     8     2
1533:       |     |     |
1534:      12--0-13--1--14
1535:     */
1536:     static Obj<ALE::Mesh> createSquareBoundary(const MPI_Comm comm, const double lower[], const double upper[], const int edges[], const int debug = 0) {
1537:       Obj<Mesh> mesh        = new Mesh(comm, 1, debug);
1538:       int       numVertices = (edges[0]+1)*(edges[1]+1);
1539:       int       numEdges    = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
1540:       double   *coords      = new double[numVertices*2];
1541:       const Obj<Mesh::sieve_type> sieve    = new Mesh::sieve_type(mesh->comm(), mesh->debug());
1542:       Mesh::point_type           *vertices = new Mesh::point_type[numVertices];
1543:       int                         order    = 0;

1545:       mesh->setSieve(sieve);
1546:       const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
1547:       if (mesh->commRank() == 0) {
1548:         /* Create sieve and ordering */
1549:         for(int v = numEdges; v < numEdges+numVertices; v++) {
1550:           vertices[v-numEdges] = Mesh::point_type(v);
1551:         }
1552:         for(int vy = 0; vy <= edges[1]; vy++) {
1553:           for(int ex = 0; ex < edges[0]; ex++) {
1554:             Mesh::point_type edge(vy*edges[0] + ex);
1555:             int vertex = vy*(edges[0]+1) + ex;

1557:             sieve->addArrow(vertices[vertex+0], edge, order++);
1558:             sieve->addArrow(vertices[vertex+1], edge, order++);
1559:             if ((vy == 0) || (vy == edges[1])) {
1560:               mesh->setValue(markers, edge, 1);
1561:               mesh->setValue(markers, vertices[vertex], 1);
1562:               if (ex == edges[0]-1) {
1563:                 mesh->setValue(markers, vertices[vertex+1], 1);
1564:               }
1565:             }
1566:           }
1567:         }
1568:         for(int vx = 0; vx <= edges[0]; vx++) {
1569:           for(int ey = 0; ey < edges[1]; ey++) {
1570:             Mesh::point_type edge(vx*edges[1] + ey + edges[0]*(edges[1]+1));
1571:             int vertex = ey*(edges[0]+1) + vx;

1573:             sieve->addArrow(vertices[vertex],            edge, order++);
1574:             sieve->addArrow(vertices[vertex+edges[0]+1], edge, order++);
1575:             if ((vx == 0) || (vx == edges[0])) {
1576:               mesh->setValue(markers, edge, 1);
1577:               mesh->setValue(markers, vertices[vertex], 1);
1578:               if (ey == edges[1]-1) {
1579:                 mesh->setValue(markers, vertices[vertex+edges[0]+1], 1);
1580:               }
1581:             }
1582:           }
1583:         }
1584:       }
1585:       mesh->stratify();
1586:       for(int vy = 0; vy <= edges[1]; ++vy) {
1587:         for(int vx = 0; vx <= edges[0]; ++vx) {
1588:           coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
1589:           coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
1590:         }
1591:       }
1592:       delete [] vertices;
1593:       ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
1594:       return mesh;
1595:     }
1598:     /*
1599:       Simple cubic boundary:

1601:      30----31-----32
1602:       |     |     |
1603:       |  3  |  2  |
1604:       |     |     |
1605:      27----28-----29
1606:       |     |     |
1607:       |  0  |  1  |
1608:       |     |     |
1609:      24----25-----26
1610:     */
1611:     static Obj<Mesh> createCubeBoundary(const MPI_Comm comm, const double lower[], const double upper[], const int faces[], const int debug = 0) {
1612:       Obj<Mesh> mesh        = new Mesh(comm, 2, debug);
1613:       int       numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
1614:       int       numFaces    = 6;
1615:       double   *coords      = new double[numVertices*3];
1616:       const Obj<Mesh::sieve_type> sieve    = new Mesh::sieve_type(mesh->comm(), mesh->debug());
1617:       Mesh::point_type           *vertices = new Mesh::point_type[numVertices];
1618:       int                         order    = 0;

1620:       mesh->setSieve(sieve);
1621:       const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
1622:       if (mesh->commRank() == 0) {
1623:         /* Create sieve and ordering */
1624:         for(int v = numFaces; v < numFaces+numVertices; v++) {
1625:           vertices[v-numFaces] = Mesh::point_type(v);
1626:           mesh->setValue(markers, vertices[v-numFaces], 1);
1627:         }
1628:         {
1629:           // Side 0 (Front)
1630:           Mesh::point_type face(0);
1631:           sieve->addArrow(vertices[0], face, order++);
1632:           sieve->addArrow(vertices[1], face, order++);
1633:           sieve->addArrow(vertices[2], face, order++);
1634:           sieve->addArrow(vertices[3], face, order++);
1635:           mesh->setValue(markers, face, 1);
1636:         }
1637:         {
1638:           // Side 1 (Back)
1639:           Mesh::point_type face(1);
1640:           sieve->addArrow(vertices[5], face, order++);
1641:           sieve->addArrow(vertices[4], face, order++);
1642:           sieve->addArrow(vertices[7], face, order++);
1643:           sieve->addArrow(vertices[6], face, order++);
1644:           mesh->setValue(markers, face, 1);
1645:         }
1646:         {
1647:           // Side 2 (Bottom)
1648:           Mesh::point_type face(2);
1649:           sieve->addArrow(vertices[4], face, order++);
1650:           sieve->addArrow(vertices[5], face, order++);
1651:           sieve->addArrow(vertices[1], face, order++);
1652:           sieve->addArrow(vertices[0], face, order++);
1653:           mesh->setValue(markers, face, 1);
1654:         }
1655:         {
1656:           // Side 3 (Top)
1657:           Mesh::point_type face(3);
1658:           sieve->addArrow(vertices[3], face, order++);
1659:           sieve->addArrow(vertices[2], face, order++);
1660:           sieve->addArrow(vertices[6], face, order++);
1661:           sieve->addArrow(vertices[7], face, order++);
1662:           mesh->setValue(markers, face, 1);
1663:         }
1664:         {
1665:           // Side 4 (Left)
1666:           Mesh::point_type face(4);
1667:           sieve->addArrow(vertices[4], face, order++);
1668:           sieve->addArrow(vertices[0], face, order++);
1669:           sieve->addArrow(vertices[3], face, order++);
1670:           sieve->addArrow(vertices[7], face, order++);
1671:           mesh->setValue(markers, face, 1);
1672:         }
1673:         {
1674:           // Side 5 (Right)
1675:           Mesh::point_type face(5);
1676:           sieve->addArrow(vertices[1], face, order++);
1677:           sieve->addArrow(vertices[5], face, order++);
1678:           sieve->addArrow(vertices[6], face, order++);
1679:           sieve->addArrow(vertices[2], face, order++);
1680:           mesh->setValue(markers, face, 1);
1681:         }
1682:       }
1683:       mesh->stratify();
1684: #if 0
1685:       for(int vz = 0; vz <= edges[2]; ++vz) {
1686:         for(int vy = 0; vy <= edges[1]; ++vy) {
1687:           for(int vx = 0; vx <= edges[0]; ++vx) {
1688:             coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
1689:             coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
1690:             coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
1691:           }
1692:         }
1693:       }
1694: #else
1695:       coords[0*3+0] = lower[0];
1696:       coords[0*3+1] = lower[1];
1697:       coords[0*3+2] = upper[2];
1698:       coords[1*3+0] = upper[0];
1699:       coords[1*3+1] = lower[1];
1700:       coords[1*3+2] = upper[2];
1701:       coords[2*3+0] = upper[0];
1702:       coords[2*3+1] = upper[1];
1703:       coords[2*3+2] = upper[2];
1704:       coords[3*3+0] = lower[0];
1705:       coords[3*3+1] = upper[1];
1706:       coords[3*3+2] = upper[2];
1707:       coords[4*3+0] = lower[0];
1708:       coords[4*3+1] = lower[1];
1709:       coords[4*3+2] = lower[2];
1710:       coords[5*3+0] = upper[0];
1711:       coords[5*3+1] = lower[1];
1712:       coords[5*3+2] = lower[2];
1713:       coords[6*3+0] = upper[0];
1714:       coords[6*3+1] = upper[1];
1715:       coords[6*3+2] = lower[2];
1716:       coords[7*3+0] = lower[0];
1717:       coords[7*3+1] = upper[1];
1718:       coords[7*3+2] = lower[2];
1719: #endif
1720:       ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
1721:       return mesh;
1722:     }
1723:   };
1724: } // namespace ALE

1726: namespace ALECompat {
1727:   using ALE::Obj;
1728:   template<typename Topology_>
1729:   class Bundle : public ALE::ParallelObject {
1730:   public:
1731:     typedef Topology_                                      topology_type;
1732:     typedef typename topology_type::patch_type             patch_type;
1733:     typedef typename topology_type::point_type             point_type;
1734:     typedef typename topology_type::sieve_type             sieve_type;
1735:     typedef typename sieve_type::coneArray                 coneArray;
1736:     typedef ALECompat::New::Section<topology_type, double> real_section_type;
1737:     typedef ALECompat::New::Section<topology_type, int>    int_section_type;
1738:     typedef ALE::MinimalArrow<point_type, point_type>      arrow_type;
1739:     typedef ALE::UniformSection<arrow_type, int>           arrow_section_type;
1740:     typedef struct {double x, y, z;}                       split_value;
1741:     typedef ALE::pair<point_type, split_value>             pair_type;
1742:     typedef ALECompat::New::Section<topology_type, pair_type>    pair_section_type;
1743:     typedef std::map<std::string, Obj<real_section_type> > real_sections_type;
1744:     typedef std::map<std::string, Obj<int_section_type> >  int_sections_type;
1745:     typedef std::map<std::string, Obj<pair_section_type> > pair_sections_type;
1746:     typedef std::map<std::string, Obj<arrow_section_type> > arrow_sections_type;
1747:     typedef typename topology_type::send_overlap_type      send_overlap_type;
1748:     typedef typename topology_type::recv_overlap_type      recv_overlap_type;
1749:     typedef typename ALECompat::New::SectionCompletion<topology_type, point_type>::topology_type      comp_topology_type;
1750:     typedef typename ALECompat::New::OverlapValues<send_overlap_type, comp_topology_type, point_type> send_section_type;
1751:     typedef typename ALECompat::New::OverlapValues<recv_overlap_type, comp_topology_type, point_type> recv_section_type;
1752:   protected:
1753:     Obj<topology_type> _topology;
1754:     bool               _distributed;
1755:     real_sections_type _realSections;
1756:     int_sections_type  _intSections;
1757:     pair_sections_type _pairSections;
1758:     arrow_sections_type _arrowSections;
1759:   public:
1760:     Bundle(MPI_Comm comm, int debug = 0) : ALE::ParallelObject(comm, debug), _distributed(false) {
1761:       this->_topology = new topology_type(comm, debug);
1762:     };
1763:     Bundle(const Obj<topology_type>& topology) : ALE::ParallelObject(topology->comm(), topology->debug()), _topology(topology), _distributed(false) {};
1764:   public: // Accessors
1765:     bool getDistributed() const {return this->_distributed;};
1766:     void setDistributed(const bool distributed) {this->_distributed = distributed;};
1767:     const Obj<topology_type>& getTopology() const {return this->_topology;};
1768:     void setTopology(const Obj<topology_type>& topology) {this->_topology = topology;};
1769:   public:
1770:     bool hasRealSection(const std::string& name) {
1771:       return this->_realSections.find(name) != this->_realSections.end();
1772:     };
1773:     const Obj<real_section_type>& getRealSection(const std::string& name) {
1774:       if (this->_realSections.find(name) == this->_realSections.end()) {
1775:         Obj<real_section_type> section = new real_section_type(this->_topology);

1777:         if (this->_debug) {std::cout << "Creating new real section: " << name << std::endl;}
1778:         this->_realSections[name] = section;
1779:       }
1780:       return this->_realSections[name];
1781:     };
1782:     void setRealSection(const std::string& name, const Obj<real_section_type>& section) {
1783:       this->_realSections[name] = section;
1784:     };
1785:     Obj<std::set<std::string> > getRealSections() const {
1786:       Obj<std::set<std::string> > names = std::set<std::string>();

1788:       for(typename real_sections_type::const_iterator s_iter = this->_realSections.begin(); s_iter != this->_realSections.end(); ++s_iter) {
1789:         names->insert(s_iter->first);
1790:       }
1791:       return names;
1792:     };
1793:     bool hasIntSection(const std::string& name) {
1794:       return this->_intSections.find(name) != this->_intSections.end();
1795:     };
1796:     const Obj<int_section_type>& getIntSection(const std::string& name) {
1797:       if (this->_intSections.find(name) == this->_intSections.end()) {
1798:         Obj<int_section_type> section = new int_section_type(this->_topology);

1800:         if (this->_debug) {std::cout << "Creating new int section: " << name << std::endl;}
1801:         this->_intSections[name] = section;
1802:       }
1803:       return this->_intSections[name];
1804:     };
1805:     void setIntSection(const std::string& name, const Obj<int_section_type>& section) {
1806:       this->_intSections[name] = section;
1807:     };
1808:     Obj<std::set<std::string> > getIntSections() const {
1809:       Obj<std::set<std::string> > names = std::set<std::string>();

1811:       for(typename int_sections_type::const_iterator s_iter = this->_intSections.begin(); s_iter != this->_intSections.end(); ++s_iter) {
1812:         names->insert(s_iter->first);
1813:       }
1814:       return names;
1815:     };
1816:     bool hasPairSection(const std::string& name) {
1817:       return this->_pairSections.find(name) != this->_pairSections.end();
1818:     };
1819:     const Obj<pair_section_type>& getPairSection(const std::string& name) {
1820:       if (this->_pairSections.find(name) == this->_pairSections.end()) {
1821:         Obj<pair_section_type> section = new pair_section_type(this->_topology);

1823:         if (this->_debug) {std::cout << "Creating new pair section: " << name << std::endl;}
1824:         this->_pairSections[name] = section;
1825:       }
1826:       return this->_pairSections[name];
1827:     };
1828:     void setPairSection(const std::string& name, const Obj<pair_section_type>& section) {
1829:       this->_pairSections[name] = section;
1830:     };
1831:     Obj<std::set<std::string> > getPairSections() const {
1832:       Obj<std::set<std::string> > names = std::set<std::string>();

1834:       for(typename pair_sections_type::const_iterator s_iter = this->_pairSections.begin(); s_iter != this->_pairSections.end(); ++s_iter) {
1835:         names->insert(s_iter->first);
1836:       }
1837:       return names;
1838:     };
1839:     bool hasArrowSection(const std::string& name) {
1840:       return this->_arrowSections.find(name) != this->_arrowSections.end();
1841:     };
1842:     const Obj<arrow_section_type>& getArrowSection(const std::string& name) {
1843:       if (this->_arrowSections.find(name) == this->_arrowSections.end()) {
1844:         Obj<arrow_section_type> section = new arrow_section_type(this->comm(), this->debug());

1846:         if (this->_debug) {std::cout << "Creating new arrow section: " << name << std::endl;}
1847:         this->_arrowSections[name] = section;
1848:       }
1849:       return this->_arrowSections[name];
1850:     };
1851:     void setArrowSection(const std::string& name, const Obj<arrow_section_type>& section) {
1852:       this->_arrowSections[name] = section;
1853:     };
1854:     Obj<std::set<std::string> > getArrowSections() const {
1855:       Obj<std::set<std::string> > names = std::set<std::string>();

1857:       for(typename arrow_sections_type::const_iterator s_iter = this->_arrowSections.begin(); s_iter != this->_arrowSections.end(); ++s_iter) {
1858:         names->insert(s_iter->first);
1859:       }
1860:       return names;
1861:     };
1862:   public: // Printing
1863:     friend std::ostream& operator<<(std::ostream& os, const split_value& s) {
1864:       os << "(" << s.x << ", "<< s.y << ", "<< s.z << ")";
1865:       return os;
1866:     };
1867:   public: // Adapter
1868:     const Obj<sieve_type>& getSieve() {return this->_topology->getPatch(0);};
1869:     int height() {return 2;};
1870:     int depth() {return 2;};
1871:   };

1873:   class Discretization : public ALE::ParallelObject {
1874:   protected:
1875:     std::map<int,int> _dim2dof;
1876:     std::map<int,int> _dim2class;
1877:   public:
1878:     Discretization(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
1879:     ~Discretization() {};
1880:   public:
1881:     const double *getQuadraturePoints() {return NULL;};
1882:     const double *getQuadratureWeights() {return NULL;};
1883:     const double *getBasis() {return NULL;};
1884:     const double *getBasisDerivatives() {return NULL;};
1885:     int  getNumDof(const int dim) {return this->_dim2dof[dim];};
1886:     void setNumDof(const int dim, const int numDof) {this->_dim2dof[dim] = numDof;};
1887:     int  getDofClass(const int dim) {return this->_dim2class[dim];};
1888:     void setDofClass(const int dim, const int dofClass) {this->_dim2class[dim] = dofClass;};
1889:   };

1891:   class BoundaryCondition : public ALE::ParallelObject {
1892:   protected:
1893:     std::string _labelName;
1894:     double    (*_func)(const double []);
1895:   public:
1896:     BoundaryCondition(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
1897:     ~BoundaryCondition() {};
1898:   public:
1899:     const std::string& getLabelName() {return this->_labelName;};
1900:     void setLabelName(const std::string& name) {this->_labelName = name;};
1901:     void setFunction(double (*func)(const double [])) {this->_func = func;};
1902:   public:
1903:     double evaluate(const double coords[]) {return this->_func(coords);};
1904:   };

1906:   class Mesh : public Bundle<ALE::Topology<int, ALE::Sieve<int,int,int> > > {
1907:   public:
1908:     typedef int                                       point_type;
1909:     typedef ALE::Sieve<point_type,int,int>            sieve_type;
1910:     typedef ALE::Topology<int, sieve_type>            topology_type;
1911:     typedef topology_type::patch_type                 patch_type;
1912:     typedef Bundle<topology_type>                     base_type;
1913:     typedef ALECompat::New::NumberingFactory<topology_type> NumberingFactory;
1914:     typedef NumberingFactory::numbering_type          numbering_type;
1915:     typedef NumberingFactory::order_type              order_type;
1916:     typedef base_type::send_overlap_type              send_overlap_type;
1917:     typedef base_type::recv_overlap_type              recv_overlap_type;
1918:     typedef base_type::send_section_type              send_section_type;
1919:     typedef base_type::recv_section_type              recv_section_type;
1920:     typedef base_type::real_sections_type             real_sections_type;
1921:     // PCICE BC
1922:     typedef struct {double rho,u,v,p;}                bc_value_type;
1923:     typedef std::map<int, bc_value_type>              bc_values_type;
1924:     // PyLith BC
1925:     typedef ALECompat::New::Section<topology_type, ALE::pair<int,double> > foliated_section_type;
1926:   protected:
1927:     int                   _dim;
1928:     Obj<NumberingFactory> _factory;
1929:     // PCICE BC
1930:     bc_values_type        _bcValues;
1931:     // PyLith BC
1932:     Obj<foliated_section_type> _boundaries;
1933:     // Discretization
1934:     Obj<Discretization>    _discretization;
1935:     Obj<BoundaryCondition> _boundaryCondition;
1936:   public:
1937:     Mesh(MPI_Comm comm, int dim, int debug = 0) : Bundle<ALE::Topology<int, ALE::Sieve<int,int,int> > >(comm, debug), _dim(dim) {
1938:       this->_factory = NumberingFactory::singleton(debug);
1939:       this->_boundaries = NULL;
1940:       this->_discretization = new Discretization(comm, debug);
1941:       this->_boundaryCondition = new BoundaryCondition(comm, debug);
1942:     };
1943:     Mesh(const Obj<topology_type>& topology, int dim) : Bundle<ALE::Topology<int, ALE::Sieve<int,int,int> > >(topology), _dim(dim) {
1944:       this->_factory = NumberingFactory::singleton(topology->debug());
1945:       this->_boundaries = NULL;
1946:     };
1947:   public: // Accessors
1948:     int getDimension() const {return this->_dim;};
1949:     void setDimension(const int dim) {this->_dim = dim;};
1950:     const Obj<NumberingFactory>& getFactory() {return this->_factory;};
1951:     const Obj<Discretization>&    getDiscretization() {return this->_discretization;};
1952:     void setDiscretization(const Obj<Discretization>& discretization) {this->_discretization = discretization;};
1953:     const Obj<BoundaryCondition>& getBoundaryCondition() {return this->_boundaryCondition;};
1954:     void setBoundaryCondition(const Obj<BoundaryCondition>& boundaryCondition) {this->_boundaryCondition = boundaryCondition;};
1955:   public: // Mesh geometry
1956: #if 0
1957:     void computeTriangleGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1958:       const double *coords = this->restrict(coordinates, e);
1959:       const int     dim    = 2;
1960:       double        invDet;

1962:       if (v0) {
1963:         for(int d = 0; d < dim; d++) {
1964:           v0[d] = coords[d];
1965:         }
1966:       }
1967:       if (J) {
1968:         for(int d = 0; d < dim; d++) {
1969:           for(int f = 0; f < dim; f++) {
1970:             J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
1971:           }
1972:         }
1973:         detJ = J[0]*J[3] - J[1]*J[2];
1974:       }
1975:       if (invJ) {
1976:         invDet  = 1.0/detJ;
1977:         invJ[0] =  invDet*J[3];
1978:         invJ[1] = -invDet*J[1];
1979:         invJ[2] = -invDet*J[2];
1980:         invJ[3] =  invDet*J[0];
1981:       }
1982:     };
1983:     static void computeTetrahedronGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1984:       const patch_type patch  = 0;
1985:       const double    *coords = coordinates->restrict(patch, e);
1986:       const int        dim    = 3;
1987:       double           invDet;

1989:       if (v0) {
1990:         for(int d = 0; d < dim; d++) {
1991:           v0[d] = coords[d];
1992:         }
1993:       }
1994:       if (J) {
1995:         for(int d = 0; d < dim; d++) {
1996:           for(int f = 0; f < dim; f++) {
1997:             J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
1998:           }
1999:         }
2000:         detJ = J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
2001:           J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
2002:           J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
2003:       }
2004:       if (invJ) {
2005:         invDet  = 1.0/detJ;
2006:         invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
2007:         invJ[0*3+1] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
2008:         invJ[0*3+2] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
2009:         invJ[1*3+0] = invDet*(J[0*3+1]*J[2*3+2] - J[0*3+2]*J[2*3+1]);
2010:         invJ[1*3+1] = invDet*(J[0*3+2]*J[2*3+0] - J[0*3+0]*J[2*3+2]);
2011:         invJ[1*3+2] = invDet*(J[0*3+0]*J[2*3+1] - J[0*3+1]*J[2*3+0]);
2012:         invJ[2*3+0] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
2013:         invJ[2*3+1] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
2014:         invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
2015:       }
2016:     };
2017:     void computeElementGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
2018:       if (this->_dim == 2) {
2019:         computeTriangleGeometry(coordinates, e, v0, J, invJ, detJ);
2020:       } else if (this->_dim == 3) {
2021:         computeTetrahedronGeometry(coordinates, e, v0, J, invJ, detJ);
2022:       } else {
2023:         throw ALE::Exception("Unsupport dimension for element geometry computation");
2024:       }
2025:     };
2026:     double getMaxVolume() {
2027:       const topology_type::sheaf_type& patches = this->getTopology()->getPatches();
2028:       double v0[3], J[9], invJ[9], detJ, maxVolume = 0.0;

2030:       for(topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2031:         const Obj<real_section_type>&             coordinates = this->getRealSection("coordinates");
2032:         const Obj<topology_type::label_sequence>& cells       = this->getTopology()->heightStratum(p_iter->first, 0);

2034:         for(topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
2035:           this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
2036:           maxVolume = std::max(maxVolume, detJ);
2037:         }
2038:       }
2039:       return maxVolume;
2040:     };
2041:     // Find the cell in which this point lies (stupid algorithm)
2042:     point_type locatePoint_2D(const patch_type& patch, const real_section_type::value_type point[]) {
2043:       const Obj<real_section_type>&             coordinates = this->getRealSection("coordinates");
2044:       const Obj<topology_type::label_sequence>& cells       = this->getTopology()->heightStratum(patch, 0);
2045:       const int                                 embedDim    = 2;
2046:       double v0[2], J[4], invJ[4], detJ;

2048:       for(topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
2049:         this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
2050:         double xi   = invJ[0*embedDim+0]*(point[0] - v0[0]) + invJ[0*embedDim+1]*(point[1] - v0[1]);
2051:         double eta  = invJ[1*embedDim+0]*(point[0] - v0[0]) + invJ[1*embedDim+1]*(point[1] - v0[1]);

2053:         if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 1.0)) {
2054:           return *c_iter;
2055:         }
2056:       }
2057:       throw ALE::Exception("Could not locate point");
2058:     };
2059:     //   Assume a simplex and 3D
2060:     point_type locatePoint_3D(const patch_type& patch, const real_section_type::value_type point[]) {
2061:       const Obj<real_section_type>&             coordinates = this->getRealSection("coordinates");
2062:       const Obj<topology_type::label_sequence>& cells       = this->getTopology()->heightStratum(patch, 0);
2063:       const int                                 embedDim    = 3;
2064:       double v0[3], J[9], invJ[9], detJ;

2066:       for(topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
2067:         this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
2068:         double xi   = invJ[0*embedDim+0]*(point[0] - v0[0]) + invJ[0*embedDim+1]*(point[1] - v0[1]) + invJ[0*embedDim+2]*(point[2] - v0[2]);
2069:         double eta  = invJ[1*embedDim+0]*(point[0] - v0[0]) + invJ[1*embedDim+1]*(point[1] - v0[1]) + invJ[1*embedDim+2]*(point[2] - v0[2]);
2070:         double zeta = invJ[2*embedDim+0]*(point[0] - v0[0]) + invJ[2*embedDim+1]*(point[1] - v0[1]) + invJ[2*embedDim+2]*(point[2] - v0[2]);

2072:         if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 1.0)) {
2073:           return *c_iter;
2074:         }
2075:       }
2076:       throw ALE::Exception("Could not locate point");
2077:     };
2078:     point_type locatePoint(const patch_type& patch, const real_section_type::value_type point[]) {
2079:       if (this->_dim == 2) {
2080:         return locatePoint_2D(patch, point);
2081:       } else if (this->_dim == 3) {
2082:         return locatePoint_3D(patch, point);
2083:       } else {
2084:         throw ALE::Exception("No point location for mesh dimension");
2085:       }
2086:     };
2087: #endif
2088:     // Only works in 2D
2089:     int orientation(const patch_type& patch, const point_type& cell) {
2090:       const Obj<real_section_type>&     coordinates = this->getRealSection("coordinates");
2091:       const Obj<topology_type>&         topology    = this->getTopology();
2092:       const Obj<sieve_type>&            sieve       = topology->getPatch(patch);
2093:       const Obj<sieve_type::coneArray>& cone        = sieve->nCone(cell, topology->depth());
2094:       sieve_type::coneArray::iterator   cBegin      = cone->begin();
2095:       real_section_type::value_type     root[2];
2096:       real_section_type::value_type     vA[2];
2097:       real_section_type::value_type     vB[2];

2099:       const real_section_type::value_type *coords = coordinates->restrictPoint(patch, *cBegin);
2100:       root[0] = coords[0];
2101:       root[1] = coords[1];
2102:       ++cBegin;
2103:       coords = coordinates->restrictPoint(patch, *cBegin);
2104:       vA[0] = coords[0] - root[0];
2105:       vA[1] = coords[1] - root[1];
2106:       ++cBegin;
2107:       coords = coordinates->restrictPoint(patch, *cBegin);
2108:       vB[0] = coords[0] - root[0];
2109:       vB[1] = coords[1] - root[1];
2110:       double det = vA[0]*vB[1] - vA[1]*vB[0];
2111:       if (det > 0.0) return  1;
2112:       if (det < 0.0) return -1;
2113:       return 0;
2114:     };
2115:   public: // BC values for PCICE
2116:     const bc_value_type& getBCValue(const int bcFunc) {
2117:       return this->_bcValues[bcFunc];
2118:     };
2119:     void setBCValue(const int bcFunc, const bc_value_type& value) {
2120:       this->_bcValues[bcFunc] = value;
2121:     };
2122:     bc_values_type& getBCValues() {
2123:       return this->_bcValues;
2124:     };
2125:     void distributeBCValues() {
2126:       int size = this->_bcValues.size();

2128:       MPI_Bcast(&size, 1, MPI_INT, 0, this->comm());
2129:       if (this->commRank()) {
2130:         for(int bc = 0; bc < size; ++bc) {
2131:           int           funcNum;
2132:           bc_value_type funcVal;

2134:           MPI_Bcast((void *) &funcNum, 1, MPI_INT,    0, this->comm());
2135:           MPI_Bcast((void *) &funcVal, 4, MPI_DOUBLE, 0, this->comm());
2136:           this->_bcValues[funcNum] = funcVal;
2137:         }
2138:       } else {
2139:         for(bc_values_type::iterator bc_iter = this->_bcValues.begin(); bc_iter != this->_bcValues.end(); ++bc_iter) {
2140:           const int&           funcNum = bc_iter->first;
2141:           const bc_value_type& funcVal = bc_iter->second;
2142:           MPI_Bcast((void *) &funcNum, 1, MPI_INT,    0, this->comm());
2143:           MPI_Bcast((void *) &funcVal, 4, MPI_DOUBLE, 0, this->comm());
2144:         }
2145:       }
2146:     };
2147:   public: // BC values for PyLith
2148:     const Obj<foliated_section_type>& getBoundariesNew() {
2149:       if (this->_boundaries.isNull()) {
2150:         this->_boundaries = new foliated_section_type(this->getTopology());
2151:       }
2152:       return this->_boundaries;
2153:     };
2154:   public: // Discretization
2155:     void setupField(const Obj<real_section_type>& s, const bool postponeGhosts = false) {
2156:       const std::string& name  = this->_boundaryCondition->getLabelName();
2157:       const patch_type   patch = 0;

2159:       for(int d = 0; d <= this->_dim; ++d) {
2160:         s->setFiberDimensionByDepth(patch, d, this->_discretization->getNumDof(d));
2161:       }
2162:       if (!name.empty()) {
2163:         const Obj<topology_type::label_sequence>& boundary = this->_topology->getLabelStratum(patch, name, 1);

2165:         for(topology_type::label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
2166:           s->setFiberDimension(patch, *e_iter, -this->_discretization->getNumDof(this->_topology->depth(patch, *e_iter)));
2167:         }
2168:       }
2169:       s->allocate(postponeGhosts);
2170:       if (!name.empty()) {
2171:         const Obj<real_section_type>&             coordinates = this->getRealSection("coordinates");
2172:         const Obj<topology_type::label_sequence>& boundary    = this->_topology->getLabelStratum(patch, name, 1);

2174:         for(topology_type::label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
2175:           const real_section_type::value_type *coords = coordinates->restrictPoint(patch, *e_iter);
2176:           const PetscScalar                    value  = this->_boundaryCondition->evaluate(coords);

2178:           s->updatePointBC(patch, *e_iter, &value);
2179:         }
2180:       }
2181:     };
2182:   public:
2183:     void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
2184:       if (comm == MPI_COMM_NULL) {
2185:         comm = this->comm();
2186:       }
2187:       if (name == "") {
2188:         PetscPrintf(comm, "viewing a Mesh\n");
2189:       } else {
2190:         PetscPrintf(comm, "viewing Mesh '%s'\n", name.c_str());
2191:       }
2192:       this->getTopology()->view("mesh topology", comm);
2193:       Obj<std::set<std::string> > sections = this->getRealSections();

2195:       for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
2196:         this->getRealSection(*name)->view(*name);
2197:       }
2198:       sections = this->getIntSections();
2199:       for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
2200:         this->getIntSection(*name)->view(*name);
2201:       }
2202:       sections = this->getPairSections();
2203:       for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
2204:         this->getPairSection(*name)->view(*name);
2205:       }
2206:     };
2207:     template<typename value_type>
2208:     static std::string printMatrix(const std::string& name, const int rows, const int cols, const value_type matrix[], const int rank = -1)
2209:     {
2210:       ostringstream output;
2211:       ostringstream rankStr;

2213:       if (rank >= 0) {
2214:         rankStr << "[" << rank << "]";
2215:       }
2216:       output << rankStr.str() << name << " = " << std::endl;
2217:       for(int r = 0; r < rows; r++) {
2218:         if (r == 0) {
2219:           output << rankStr.str() << " /";
2220:         } else if (r == rows-1) {
2221:           output << rankStr.str() << " \\";
2222:         } else {
2223:           output << rankStr.str() << " |";
2224:         }
2225:         for(int c = 0; c < cols; c++) {
2226:           output << " " << matrix[r*cols+c];
2227:         }
2228:         if (r == 0) {
2229:           output << " \\" << std::endl;
2230:         } else if (r == rows-1) {
2231:           output << " /" << std::endl;
2232:         } else {
2233:           output << " |" << std::endl;
2234:         }
2235:       }
2236:       return output.str();
2237:     };
2238:   };
2239: } // namespace ALECompat
2240: #endif