Actual source code: Mesh.hh
1: #ifndef included_ALE_Mesh_hh
2: #define included_ALE_Mesh_hh
4: #include <valarray>
6: #ifndef included_ALE_Numbering_hh
7: #include <Numbering.hh>
8: #endif
10: #ifndef included_ALE_Field_hh
11: #include <Field.hh>
12: #endif
14: #ifndef included_ALE_SieveBuilder_hh
15: #include <SieveBuilder.hh>
16: #endif
18: #ifndef included_ALE_LabelSifter_hh
19: #include <LabelSifter.hh>
20: #endif
22: namespace ALE {
23: class indexSet : public std::valarray<int> {
24: public:
25: inline bool
26: operator<(const indexSet& __x) {
27: if (__x.size() != this->size()) return __x.size() < this->size();
28: for(unsigned int i = 0; i < __x.size(); ++i) {
29: if (__x[i] == (*this)[i]) continue;
30: return __x[i] < (*this)[i];
31: }
32: return false;
33: }
34: };
35: inline bool
36: operator<(const indexSet& __x, const indexSet& __y) {
37: if (__x.size() != __y.size()) return __x.size() < __y.size();
38: for(unsigned int i = 0; i < __x.size(); ++i) {
39: if (__x[i] == __y[i]) continue;
40: return __x[i] < __y[i];
41: }
42: return false;
43: };
44: inline bool
45: operator<=(const indexSet& __x, const indexSet& __y) {
46: if (__x.size() != __y.size()) return __x.size() < __y.size();
47: for(unsigned int i = 0; i < __x.size(); ++i) {
48: if (__x[i] == __y[i]) continue;
49: return __x[i] < __y[i];
50: }
51: return true;
52: };
53: inline bool
54: operator==(const indexSet& __x, const indexSet& __y) {
55: if (__x.size() != __y.size()) return false;
56: for(unsigned int i = 0; i < __x.size(); ++i) {
57: if (__x[i] != __y[i]) return false;
58: }
59: return true;
60: };
61: inline bool
62: operator!=(const indexSet& __x, const indexSet& __y) {
63: if (__x.size() != __y.size()) return true;
64: for(unsigned int i = 0; i < __x.size(); ++i) {
65: if (__x[i] != __y[i]) return true;
66: }
67: return false;
68: };
70: template<typename Sieve_,
71: typename RealSection_ = Section<typename Sieve_::point_type, double>,
72: typename IntSection_ = Section<typename Sieve_::point_type, int>,
73: typename ArrowSection_ = UniformSection<MinimalArrow<typename Sieve_::point_type, typename Sieve_::point_type>, int> >
74: class Bundle : public ALE::ParallelObject {
75: public:
76: typedef Sieve_ sieve_type;
77: typedef RealSection_ real_section_type;
78: typedef IntSection_ int_section_type;
79: typedef ArrowSection_ arrow_section_type;
80: typedef Bundle<Sieve_,RealSection_,IntSection_,ArrowSection_> this_type;
81: typedef typename sieve_type::point_type point_type;
82: typedef malloc_allocator<point_type> alloc_type;
83: #define NEW_LABEL
84: #ifdef NEW_LABEL
85: typedef typename ALE::LabelSifter<int, point_type> label_type;
86: #else
87: typedef typename ALE::Sifter<int, point_type, int> label_type;
88: #endif
89: typedef typename std::map<const std::string, Obj<label_type> > labels_type;
90: typedef typename label_type::supportSequence label_sequence;
91: typedef std::map<std::string, Obj<arrow_section_type> > arrow_sections_type;
92: typedef std::map<std::string, Obj<real_section_type> > real_sections_type;
93: typedef std::map<std::string, Obj<int_section_type> > int_sections_type;
94: typedef ALE::Point index_type;
95: typedef std::pair<index_type, int> oIndex_type;
96: typedef std::vector<oIndex_type> oIndexArray;
97: typedef std::pair<int *, int> indices_type;
98: typedef NumberingFactory<this_type> NumberingFactory;
99: typedef typename NumberingFactory::numbering_type numbering_type;
100: typedef typename NumberingFactory::order_type order_type;
101: typedef typename ALE::Sifter<int,point_type,point_type> send_overlap_type;
102: typedef typename ALE::Sifter<point_type,int,point_type> recv_overlap_type;
103: typedef typename ALE::SieveAlg<this_type> sieve_alg_type;
104: typedef typename sieve_alg_type::coneArray coneArray;
105: typedef typename sieve_alg_type::orientedConeArray oConeArray;
106: typedef typename sieve_alg_type::supportArray supportArray;
107: protected:
108: Obj<sieve_type> _sieve;
109: labels_type _labels;
110: int _maxHeight;
111: int _maxDepth;
112: arrow_sections_type _arrowSections;
113: real_sections_type _realSections;
114: int_sections_type _intSections;
115: Obj<oIndexArray> _indexArray;
116: Obj<NumberingFactory> _factory;
117: bool _calculatedOverlap;
118: Obj<send_overlap_type> _sendOverlap;
119: Obj<recv_overlap_type> _recvOverlap;
120: Obj<send_overlap_type> _distSendOverlap;
121: Obj<recv_overlap_type> _distRecvOverlap;
122: // Work space
123: Obj<std::set<point_type> > _modifiedPoints;
124: public:
125: Bundle(MPI_Comm comm, int debug = 0) : ALE::ParallelObject(comm, debug), _maxHeight(-1), _maxDepth(-1) {
126: this->_indexArray = new oIndexArray();
127: this->_modifiedPoints = new std::set<point_type>();
128: this->_factory = NumberingFactory::singleton(this->comm(), this->debug());
129: this->_calculatedOverlap = false;
130: this->_sendOverlap = new send_overlap_type(comm, debug);
131: this->_recvOverlap = new recv_overlap_type(comm, debug);
132: };
133: Bundle(const Obj<sieve_type>& sieve) : ALE::ParallelObject(sieve->comm(), sieve->debug()), _sieve(sieve), _maxHeight(-1), _maxDepth(-1) {
134: this->_indexArray = new oIndexArray();
135: this->_modifiedPoints = new std::set<point_type>();
136: this->_factory = NumberingFactory::singleton(this->comm(), this->debug());
137: this->_calculatedOverlap = false;
138: this->_sendOverlap = new send_overlap_type(comm, debug);
139: this->_recvOverlap = new recv_overlap_type(comm, debug);
140: };
141: virtual ~Bundle() {};
142: public: // Verifiers
143: bool hasLabel(const std::string& name) {
144: if (this->_labels.find(name) != this->_labels.end()) {
145: return true;
146: }
147: return false;
148: };
149: void checkLabel(const std::string& name) {
150: if (!this->hasLabel(name)) {
151: ostringstream msg;
152: msg << "Invalid label name: " << name << std::endl;
153: throw ALE::Exception(msg.str().c_str());
154: }
155: };
156: public: // Accessors
157: const Obj<sieve_type>& getSieve() const {return this->_sieve;};
158: void setSieve(const Obj<sieve_type>& sieve) {this->_sieve = sieve;};
159: bool hasArrowSection(const std::string& name) const {
160: return this->_arrowSections.find(name) != this->_arrowSections.end();
161: };
162: const Obj<arrow_section_type>& getArrowSection(const std::string& name) {
163: if (!this->hasArrowSection(name)) {
164: Obj<arrow_section_type> section = new arrow_section_type(this->comm(), this->debug());
166: section->setName(name);
167: if (this->_debug) {std::cout << "Creating new arrow section: " << name << std::endl;}
168: this->_arrowSections[name] = section;
169: }
170: return this->_arrowSections[name];
171: };
172: void setArrowSection(const std::string& name, const Obj<arrow_section_type>& section) {
173: this->_arrowSections[name] = section;
174: };
175: Obj<std::set<std::string> > getArrowSections() const {
176: Obj<std::set<std::string> > names = std::set<std::string>();
178: for(typename arrow_sections_type::const_iterator s_iter = this->_arrowSections.begin(); s_iter != this->_arrowSections.end(); ++s_iter) {
179: names->insert(s_iter->first);
180: }
181: return names;
182: };
183: bool hasRealSection(const std::string& name) const {
184: return this->_realSections.find(name) != this->_realSections.end();
185: };
186: const Obj<real_section_type>& getRealSection(const std::string& name) {
187: if (!this->hasRealSection(name)) {
188: Obj<real_section_type> section = new real_section_type(this->comm(), this->debug());
190: section->setName(name);
191: if (this->_debug) {std::cout << "Creating new real section: " << name << std::endl;}
192: this->_realSections[name] = section;
193: }
194: return this->_realSections[name];
195: };
196: void setRealSection(const std::string& name, const Obj<real_section_type>& section) {
197: this->_realSections[name] = section;
198: };
199: Obj<std::set<std::string> > getRealSections() const {
200: Obj<std::set<std::string> > names = std::set<std::string>();
202: for(typename real_sections_type::const_iterator s_iter = this->_realSections.begin(); s_iter != this->_realSections.end(); ++s_iter) {
203: names->insert(s_iter->first);
204: }
205: return names;
206: };
207: bool hasIntSection(const std::string& name) const {
208: return this->_intSections.find(name) != this->_intSections.end();
209: };
210: const Obj<int_section_type>& getIntSection(const std::string& name) {
211: if (!this->hasIntSection(name)) {
212: Obj<int_section_type> section = new int_section_type(this->comm(), this->debug());
214: section->setName(name);
215: if (this->_debug) {std::cout << "Creating new int section: " << name << std::endl;}
216: this->_intSections[name] = section;
217: }
218: return this->_intSections[name];
219: };
220: void setIntSection(const std::string& name, const Obj<int_section_type>& section) {
221: this->_intSections[name] = section;
222: };
223: Obj<std::set<std::string> > getIntSections() const {
224: Obj<std::set<std::string> > names = std::set<std::string>();
226: for(typename int_sections_type::const_iterator s_iter = this->_intSections.begin(); s_iter != this->_intSections.end(); ++s_iter) {
227: names->insert(s_iter->first);
228: }
229: return names;
230: };
231: const Obj<NumberingFactory>& getFactory() const {return this->_factory;};
232: const Obj<send_overlap_type>& getSendOverlap() const {return this->_sendOverlap;};
233: void setSendOverlap(const Obj<send_overlap_type>& overlap) {this->_sendOverlap = overlap;};
234: const Obj<recv_overlap_type>& getRecvOverlap() const {return this->_recvOverlap;};
235: void setRecvOverlap(const Obj<recv_overlap_type>& overlap) {this->_recvOverlap = overlap;};
236: const Obj<send_overlap_type>& getDistSendOverlap() const {return this->_distSendOverlap;};
237: void setDistSendOverlap(const Obj<send_overlap_type>& overlap) {this->_distSendOverlap = overlap;};
238: const Obj<recv_overlap_type>& getDistRecvOverlap() const {return this->_distRecvOverlap;};
239: void setDistRecvOverlap(const Obj<recv_overlap_type>& overlap) {this->_distRecvOverlap = overlap;};
240: public: // Labels
241: int getValue (const Obj<label_type>& label, const point_type& point, const int defValue = 0) {
242: const Obj<typename label_type::coneSequence>& cone = label->cone(point);
244: if (cone->size() == 0) return defValue;
245: return *cone->begin();
246: };
247: void setValue(const Obj<label_type>& label, const point_type& point, const int value) {
248: label->setCone(value, point);
249: };
250: template<typename InputPoints>
251: int getMaxValue (const Obj<label_type>& label, const Obj<InputPoints>& points, const int defValue = 0) {
252: int maxValue = defValue;
254: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
255: maxValue = std::max(maxValue, this->getValue(label, *p_iter, defValue));
256: }
257: return maxValue;
258: };
259: const Obj<label_type>& createLabel(const std::string& name) {
260: this->_labels[name] = new label_type(this->comm(), this->debug());
261: return this->_labels[name];
262: };
263: const Obj<label_type>& getLabel(const std::string& name) {
264: this->checkLabel(name);
265: return this->_labels[name];
266: };
267: void setLabel(const std::string& name, const Obj<label_type>& label) {
268: this->_labels[name] = label;
269: };
270: const labels_type& getLabels() {
271: return this->_labels;
272: };
273: virtual const Obj<label_sequence>& getLabelStratum(const std::string& name, int value) {
274: this->checkLabel(name);
275: return this->_labels[name]->support(value);
276: };
277: public: // Stratification
278: template<class InputPoints>
279: void computeHeight(const Obj<label_type>& height, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxHeight) {
280: this->_modifiedPoints->clear();
282: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
283: // Compute the max height of the points in the support of p, and add 1
284: int h0 = this->getValue(height, *p_iter, -1);
285: int h1 = this->getMaxValue(height, sieve->support(*p_iter), -1) + 1;
287: if(h1 != h0) {
288: this->setValue(height, *p_iter, h1);
289: if (h1 > maxHeight) maxHeight = h1;
290: this->_modifiedPoints->insert(*p_iter);
291: }
292: }
293: // FIX: We would like to avoid the copy here with cone()
294: if(this->_modifiedPoints->size() > 0) {
295: this->computeHeight(height, sieve, sieve->cone(this->_modifiedPoints), maxHeight);
296: }
297: };
298: void computeHeights() {
299: const Obj<label_type>& label = this->createLabel(std::string("height"));
301: this->_maxHeight = -1;
302: this->computeHeight(label, this->_sieve, this->_sieve->leaves(), this->_maxHeight);
303: };
304: virtual int height() const {return this->_maxHeight;};
305: virtual int height(const point_type& point) {
306: return this->getValue(this->_labels["height"], point, -1);
307: };
308: virtual const Obj<label_sequence>& heightStratum(int height) {
309: return this->getLabelStratum("height", height);
310: };
311: void setHeight(const Obj<label_type>& label) {
312: this->_labels["height"] = label;
313: const Obj<typename label_type::traits::capSequence> cap = label->cap();
315: for(typename label_type::traits::capSequence::iterator c_iter = cap->begin(); c_iter != cap->end(); ++c_iter) {
316: this->_maxHeight = std::max(this->_maxHeight, *c_iter);
317: }
318: };
319: template<class InputPoints>
320: void computeDepth(const Obj<label_type>& depth, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxDepth) {
321: this->_modifiedPoints->clear();
323: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
324: // Compute the max depth of the points in the cone of p, and add 1
325: int d0 = this->getValue(depth, *p_iter, -1);
326: int d1 = this->getMaxValue(depth, sieve->cone(*p_iter), -1) + 1;
328: if(d1 != d0) {
329: this->setValue(depth, *p_iter, d1);
330: if (d1 > maxDepth) maxDepth = d1;
331: this->_modifiedPoints->insert(*p_iter);
332: }
333: }
334: // FIX: We would like to avoid the copy here with support()
335: if(this->_modifiedPoints->size() > 0) {
336: this->computeDepth(depth, sieve, sieve->support(this->_modifiedPoints), maxDepth);
337: }
338: };
339: void computeDepths() {
340: const Obj<label_type>& label = this->createLabel(std::string("depth"));
342: this->_maxDepth = -1;
343: this->computeDepth(label, this->_sieve, this->_sieve->roots(), this->_maxDepth);
344: };
345: virtual int depth() const {return this->_maxDepth;};
346: virtual int depth(const point_type& point) {
347: return this->getValue(this->_labels["depth"], point, -1);
348: };
349: virtual const Obj<label_sequence>& depthStratum(int depth) {
350: return this->getLabelStratum("depth", depth);
351: };
354: virtual void stratify() {
355: ALE_LOG_EVENT_BEGIN;
356: this->computeHeights();
357: this->computeDepths();
358: ALE_LOG_EVENT_END;
359: };
360: public: // Size traversal
361: template<typename Section_>
362: int size(const Obj<Section_>& section, const point_type& p) {
363: const typename Section_::chart_type& chart = section->getChart();
364: int size = 0;
366: if (this->height() < 2) {
367: const Obj<typename sieve_type::coneSequence>& cone = this->_sieve->cone(p);
368: typename sieve_type::coneSequence::iterator end = cone->end();
370: if (chart.count(p)) {
371: size += section->getConstrainedFiberDimension(p);
372: }
373: for(typename sieve_type::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
374: if (chart.count(*c_iter)) {
375: size += section->getConstrainedFiberDimension(*c_iter);
376: }
377: }
378: } else {
379: const Obj<coneArray> closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), p);
380: typename coneArray::iterator end = closure->end();
382: for(typename coneArray::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
383: if (chart.count(*c_iter)) {
384: size += section->getConstrainedFiberDimension(*c_iter);
385: }
386: }
387: }
388: return size;
389: };
390: template<typename Section_>
391: int sizeWithBC(const Obj<Section_>& section, const point_type& p) {
392: const typename Section_::chart_type& chart = section->getChart();
393: int size = 0;
395: if (this->height() < 2) {
396: const Obj<typename sieve_type::coneSequence>& cone = this->_sieve->cone(p);
397: typename sieve_type::coneSequence::iterator end = cone->end();
399: if (chart.count(p)) {
400: size += section->getFiberDimension(p);
401: }
402: for(typename sieve_type::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
403: if (chart.count(*c_iter)) {
404: size += section->getFiberDimension(*c_iter);
405: }
406: }
407: } else {
408: const Obj<coneArray> closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), p);
409: typename coneArray::iterator end = closure->end();
411: for(typename coneArray::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
412: if (chart.count(*c_iter)) {
413: size += section->getFiberDimension(*c_iter);
414: }
415: }
416: }
417: return size;
418: };
419: protected:
420: int *getIndexArray(const int size) {
421: static int *array = NULL;
422: static int maxSize = 0;
424: if (size > maxSize) {
425: maxSize = size;
426: if (array) delete [] array;
427: array = new int[maxSize];
428: };
429: return array;
430: };
431: public: // Index traversal
432: void expandInterval(const index_type& interval, PetscInt indices[], PetscInt *indx) {
433: const int end = interval.prefix + interval.index;
435: for(int i = interval.index; i < end; ++i) {
436: indices[(*indx)++] = i;
437: }
438: };
439: void expandInterval(const index_type& interval, const int orientation, PetscInt indices[], PetscInt *indx) {
440: if (orientation >= 0) {
441: for(int i = 0; i < interval.prefix; ++i) {
442: indices[(*indx)++] = interval.index + i;
443: }
444: } else {
445: for(int i = interval.prefix-1; i >= 0; --i) {
446: indices[(*indx)++] = interval.index + i;
447: }
448: }
449: for(int i = 0; i < -interval.prefix; ++i) {
450: indices[(*indx)++] = -1;
451: }
452: };
453: void expandIntervals(Obj<oIndexArray> intervals, PetscInt *indices) {
454: int k = 0;
456: for(typename oIndexArray::iterator i_iter = intervals->begin(); i_iter != intervals->end(); ++i_iter) {
457: this->expandInterval(i_iter->first, i_iter->second, indices, &k);
458: }
459: }
460: template<typename Section_>
461: const indices_type getIndicesRaw(const Obj<Section_>& section, const point_type& p) {
462: int *indexArray = NULL;
463: int size = 0;
465: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
466: typename oConeArray::iterator begin = closure->begin();
467: typename oConeArray::iterator end = closure->end();
469: for(typename oConeArray::iterator p_iter = begin; p_iter != end; ++p_iter) {
470: size += section->getFiberDimension(p_iter->first);
471: }
472: indexArray = this->getIndexArray(size);
473: int k = 0;
474: for(typename oConeArray::iterator p_iter = begin; p_iter != end; ++p_iter) {
475: section->getIndicesRaw(p_iter->first, section->getIndex(p_iter->first), indexArray, &k, p_iter->second);
476: }
477: return indices_type(indexArray, size);
478: };
479: template<typename Section_>
480: const indices_type getIndices(const Obj<Section_>& section, const point_type& p, const int level = -1) {
481: int *indexArray = NULL;
482: int size = 0;
484: if (level == 0) {
485: size += section->getFiberDimension(p);
486: indexArray = this->getIndexArray(size);
487: int k = 0;
489: section->getIndices(p, indexArray, &k);
490: } else if ((level == 1) || (this->height() == 1)) {
491: const Obj<typename sieve_type::coneSequence>& cone = this->_sieve->cone(p);
492: typename sieve_type::coneSequence::iterator end = cone->end();
494: size += section->getFiberDimension(p);
495: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
496: size += section->getFiberDimension(*p_iter);
497: }
498: indexArray = this->getIndexArray(size);
499: int k = 0;
501: section->getIndices(p, indexArray, &k);
502: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
503: section->getIndices(*p_iter, indexArray, &k);
504: }
505: } else if (level == -1) {
506: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
507: typename oConeArray::iterator begin = closure->begin();
508: typename oConeArray::iterator end = closure->end();
510: for(typename oConeArray::iterator p_iter = begin; p_iter != end; ++p_iter) {
511: size += section->getFiberDimension(p_iter->first);
512: }
513: indexArray = this->getIndexArray(size);
514: int k = 0;
515: for(typename oConeArray::iterator p_iter = begin; p_iter != end; ++p_iter) {
516: section->getIndices(p_iter->first, indexArray, &k, p_iter->second);
517: }
518: } else {
519: throw ALE::Exception("Bundle has not yet implemented getIndices() for an arbitrary level");
520: }
521: if (this->debug()) {
522: for(int i = 0; i < size; ++i) {
523: printf("[%d]index %d: %d\n", this->commRank(), i, indexArray[i]);
524: }
525: }
526: return indices_type(indexArray, size);
527: };
528: template<typename Section_, typename Numbering>
529: const indices_type getIndices(const Obj<Section_>& section, const point_type& p, const Obj<Numbering>& numbering, const int level = -1) {
530: int *indexArray = NULL;
531: int size = 0;
533: if (level == 0) {
534: size += section->getFiberDimension(p);
535: indexArray = this->getIndexArray(size);
536: int k = 0;
538: section->getIndices(p, numbering, indexArray, &k);
539: } else if ((level == 1) || (this->height() == 1)) {
540: const Obj<typename sieve_type::coneSequence>& cone = this->_sieve->cone(p);
541: typename sieve_type::coneSequence::iterator end = cone->end();
543: size += section->getFiberDimension(p);
544: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
545: size += section->getFiberDimension(*p_iter);
546: }
547: indexArray = this->getIndexArray(size);
548: int k = 0;
550: section->getIndices(p, numbering, indexArray, &k);
551: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
552: section->getIndices(*p_iter, numbering, indexArray, &k);
553: }
554: } else if (level == -1) {
555: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
556: typename oConeArray::iterator end = closure->end();
558: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
559: size += section->getFiberDimension(p_iter->first);
560: }
561: indexArray = this->getIndexArray(size);
562: int k = 0;
563: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
564: section->getIndices(p_iter->first, numbering, indexArray, &k, p_iter->second);
565: }
566: } else {
567: throw ALE::Exception("Bundle has not yet implemented getIndices() for an arbitrary level");
568: }
569: return indices_type(indexArray, size);
570: };
571: public: // Retrieval traversal
572: // Return the values for the closure of this point
573: // use a smart pointer?
574: template<typename Section_>
575: const typename Section_::value_type *restrict(const Obj<Section_>& section, const point_type& p) {
576: const int size = this->sizeWithBC(section, p);
577: return this->restrict(section, p, section->getRawArray(size), size);
578: };
579: template<typename Section_>
580: const typename Section_::value_type *restrict(const Obj<Section_>& section, const point_type& p, typename Section_::value_type *values, const int valuesSize) {
581: const int size = this->sizeWithBC(section, p);
582: int j = -1;
583: if (valuesSize < size) throw ALE::Exception("Input array too small");
585: // We could actually ask for the height of the individual point
586: if (this->height() < 2) {
587: const int& dim = section->getFiberDimension(p);
588: const typename Section_::value_type *array = section->restrictPoint(p);
590: for(int i = 0; i < dim; ++i) {
591: values[++j] = array[i];
592: }
593: const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
594: typename sieve_type::coneSequence::iterator end = cone->end();
596: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
597: const int& dim = section->getFiberDimension(*p_iter);
599: array = section->restrictPoint(*p_iter);
600: for(int i = 0; i < dim; ++i) {
601: values[++j] = array[i];
602: }
603: }
604: } else {
605: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
606: typename oConeArray::iterator end = closure->end();
608: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
609: const typename Section_::value_type *array = section->restrictPoint(p_iter->first);
610: const int& dim = section->getFiberDimension(p_iter->first);
612: if (p_iter->second >= 0) {
613: for(int i = 0; i < dim; ++i) {
614: values[++j] = array[i];
615: }
616: } else {
617: for(int i = dim-1; i >= 0; --i) {
618: values[++j] = array[i];
619: }
620: }
621: }
622: }
623: if (j != size-1) {
624: ostringstream txt;
626: txt << "Invalid restrict to point " << p << std::endl;
627: txt << " j " << j << " should be " << (size-1) << std::endl;
628: std::cout << txt.str();
629: throw ALE::Exception(txt.str().c_str());
630: }
631: return values;
632: };
633: template<typename Section_>
634: const typename Section_::value_type *restrictNew(const Obj<Section_>& section, const point_type& p) {
635: const int size = this->sizeWithBC(section, p);
636: return this->restrictNew(section, p, section->getRawArray(size), size);
637: };
638: template<typename Section_>
639: const typename Section_::value_type *restrictNew(const Obj<Section_>& section, const point_type& p, typename Section_::value_type *values, const int valuesSize) {
640: const int size = this->sizeWithBC(section, p);
641: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
642: typename oConeArray::iterator end = closure->end();
643: int j = -1;
644: if (valuesSize < size) throw ALE::Exception("Input array too small");
646: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
647: const typename Section_::value_type *array = section->restrictPoint(p_iter->first);
649: if (p_iter->second >= 0) {
650: const int& dim = section->getFiberDimension(p_iter->first);
652: for(int i = 0; i < dim; ++i) {
653: values[++j] = array[i];
654: }
655: } else {
656: int offset = 0;
658: for(int space = 0; space < section->getNumSpaces(); ++space) {
659: const int& dim = section->getFiberDimension(p_iter->first, space);
661: for(int i = dim-1; i >= 0; --i) {
662: values[++j] = array[i+offset];
663: }
664: offset += dim;
665: }
666: }
667: }
668: if (j != size-1) {
669: ostringstream txt;
671: txt << "Invalid restrict to point " << p << std::endl;
672: txt << " j " << j << " should be " << (size-1) << std::endl;
673: std::cout << txt.str();
674: throw ALE::Exception(txt.str().c_str());
675: }
676: return values;
677: };
678: template<typename Section_>
679: void update(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
680: int j = 0;
682: if (this->height() < 2) {
683: section->updatePoint(p, &v[j]);
684: j += section->getFiberDimension(p);
685: const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
687: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
688: section->updatePoint(*p_iter, &v[j]);
689: j += section->getFiberDimension(*p_iter);
690: }
691: } else {
692: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
693: typename oConeArray::iterator end = closure->end();
695: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
696: section->updatePoint(p_iter->first, &v[j], p_iter->second);
697: j += section->getFiberDimension(p_iter->first);
698: }
699: }
700: };
701: template<typename Section_>
702: void updateAdd(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
703: int j = 0;
705: if (this->height() < 2) {
706: section->updateAddPoint(p, &v[j]);
707: j += section->getFiberDimension(p);
708: const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
710: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
711: section->updateAddPoint(*p_iter, &v[j]);
712: j += section->getFiberDimension(*p_iter);
713: }
714: } else {
715: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
716: typename oConeArray::iterator end = closure->end();
718: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
719: section->updateAddPoint(p_iter->first, &v[j], p_iter->second);
720: j += section->getFiberDimension(p_iter->first);
721: }
722: }
723: };
724: template<typename Section_>
725: void updateBC(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
726: int j = 0;
728: if (this->height() < 2) {
729: section->updatePointBC(p, &v[j]);
730: j += section->getFiberDimension(p);
731: const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
733: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
734: section->updatePointBC(*p_iter, &v[j]);
735: j += section->getFiberDimension(*p_iter);
736: }
737: } else {
738: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
739: typename oConeArray::iterator end = closure->end();
741: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
742: section->updatePointBC(p_iter->first, &v[j], p_iter->second);
743: j += section->getFiberDimension(p_iter->first);
744: }
745: }
746: };
747: template<typename Section_>
748: void updateAll(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
749: int j = 0;
751: if (this->height() < 2) {
752: section->updatePointAll(p, &v[j]);
753: j += section->getFiberDimension(p);
754: const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
756: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
757: section->updatePointAll(*p_iter, &v[j]);
758: j += section->getFiberDimension(*p_iter);
759: }
760: } else {
761: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
762: typename oConeArray::iterator end = closure->end();
764: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
765: section->updatePointAll(p_iter->first, &v[j], p_iter->second);
766: j += section->getFiberDimension(p_iter->first);
767: }
768: }
769: };
770: template<typename Section_>
771: void updateAllAdd(const Obj<Section_>& section, const point_type& p, const typename Section_::value_type v[]) {
772: int j = 0;
774: if (this->height() < 2) {
775: section->updatePointAllAdd(p, &v[j]);
776: j += section->getFiberDimension(p);
777: const Obj<typename sieve_type::coneSequence>& cone = this->getSieve()->cone(p);
779: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
780: section->updatePointAllAdd(*p_iter, &v[j]);
781: j += section->getFiberDimension(*p_iter);
782: }
783: } else {
784: const Obj<oConeArray> closure = sieve_alg_type::orientedClosure(this, this->getArrowSection("orientation"), p);
785: typename oConeArray::iterator end = closure->end();
787: for(typename oConeArray::iterator p_iter = closure->begin(); p_iter != end; ++p_iter) {
788: section->updatePointAllAdd(p_iter->first, &v[j], p_iter->second);
789: j += section->getFiberDimension(p_iter->first);
790: }
791: }
792: };
793: public: // Optimization
794: // Calculate a custom atlas for the given traversal
795: // This returns the tag value assigned to the traversal
796: template<typename Section_, typename Sequence_>
797: int calculateCustomAtlas(const Obj<Section_>& section, const Obj<Sequence_>& points) {
798: const typename Sequence_::iterator begin = points->begin();
799: const typename Sequence_::iterator end = points->end();
800: const int num = points->size();
801: int *rOffsets = new int[num+1];
802: int *rIndices;
803: int *uOffsets = new int[num+1];
804: int *uIndices;
805: int p;
807: p = 0;
808: rOffsets[p] = 0;
809: uOffsets[p] = 0;
810: for(typename Sequence_::iterator p_iter = begin; p_iter != end; ++p_iter, ++p) {
811: rOffsets[p+1] = rOffsets[p] + this->sizeWithBC(section, *p_iter);
812: uOffsets[p+1] = rOffsets[p+1];
813: //uOffsets[p+1] = uOffsets[p] + this->size(section, *p_iter);
814: }
815: rIndices = new int[rOffsets[p]];
816: uIndices = new int[uOffsets[p]];
817: p = 0;
818: for(typename Sequence_::iterator p_iter = begin; p_iter != end; ++p_iter, ++p) {
819: const indices_type rIdx = this->getIndicesRaw(section, *p_iter);
820: for(int i = 0, k = rOffsets[p]; k < rOffsets[p+1]; ++i, ++k) rIndices[k] = rIdx.first[i];
822: const indices_type uIdx = this->getIndices(section, *p_iter);
823: for(int i = 0, k = uOffsets[p]; k < uOffsets[p+1]; ++i, ++k) uIndices[k] = uIdx.first[i];
824: }
825: return section->setCustomAtlas(rOffsets, rIndices, uOffsets, uIndices);
826: };
827: template<typename Section_>
828: const typename Section_::value_type *restrict(const Obj<Section_>& section, const int tag, const int p) {
829: const int *offsets, *indices;
831: section->getCustomRestrictAtlas(tag, &offsets, &indices);
832: const int size = offsets[p+1] - offsets[p];
833: return this->restrict(section, tag, p, section->getRawArray(size), offsets, indices);
834: };
835: template<typename Section_>
836: const typename Section_::value_type *restrict(const Obj<Section_>& section, const int tag, const int p, typename Section_::value_type *values, const int valuesSize) {
837: const int *offsets, *indices;
839: section->getCustomRestrictAtlas(tag, &offsets, &indices);
840: const int size = offsets[p+1] - offsets[p];
841: if (valuesSize < size) {throw ALE::Exception("Input array too small");}
842: return this->restrict(section, tag, p, values, offsets, indices);
843: };
844: template<typename Section_>
845: const typename Section_::value_type *restrict(const Obj<Section_>& section, const int tag, const int p, typename Section_::value_type *values, const int offsets[], const int indices[]) {
846: const typename Section_::value_type *array = section->restrict();
848: const int size = offsets[p+1] - offsets[p];
849: for(int j = 0, k = offsets[p]; j < size; ++j, ++k) {
850: values[j] = array[indices[k]];
851: }
852: return values;
853: };
854: template<typename Section_>
855: void updateAdd(const Obj<Section_>& section, const int tag, const int p, const typename Section_::value_type values[]) {
856: typename Section_::value_type *array = (typename Section_::value_type *) section->restrict();
857: const int *offsets, *indices;
859: section->getCustomUpdateAtlas(tag, &offsets, &indices);
860: const int size = offsets[p+1] - offsets[p];
861: for(int j = 0, k = offsets[p]; j < size; ++j, ++k) {
862: if (indices[k] < 0) continue;
863: array[indices[k]] += values[j];
864: }
865: };
866: public: // Allocation
867: template<typename Section_>
868: void allocate(const Obj<Section_>& section, const Obj<send_overlap_type>& sendOverlap = NULL) {
869: bool doGhosts = !sendOverlap.isNull();
871: this->_factory->orderPatch(section, this->getSieve(), sendOverlap);
872: if (doGhosts) {
873: if (this->_debug > 1) {std::cout << "Ordering patch for ghosts" << std::endl;}
874: const typename Section_::chart_type& points = section->getChart();
875: typename Section_::index_type::index_type offset = 0;
877: for(typename Section_::chart_type::const_iterator point = points.begin(); point != points.end(); ++point) {
878: const typename Section_::index_type& idx = section->getIndex(*point);
880: offset = std::max(offset, idx.index + std::abs(idx.prefix));
881: }
882: this->_factory->orderPatch(section, this->getSieve(), NULL, offset);
883: if (offset != section->sizeWithBC()) throw ALE::Exception("Inconsistent array sizes in section");
884: }
885: section->allocateStorage();
886: };
887: template<typename Section_>
888: void reallocate(const Obj<Section_>& section) {
889: if (section->getNewAtlas().isNull()) return;
890: // Since copy() preserves offsets, we must reinitialize them before ordering
891: const Obj<typename Section_::atlas_type> atlas = section->getAtlas();
892: const Obj<typename Section_::atlas_type>& newAtlas = section->getNewAtlas();
893: const typename Section_::atlas_type::chart_type& chart = newAtlas->getChart();
894: const typename Section_::atlas_type::chart_type& oldChart = atlas->getChart();
895: int newSize = 0;
896: typename Section_::index_type defaultIdx(0, -1);
898: for(typename Section_::atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
899: defaultIdx.prefix = newAtlas->restrictPoint(*c_iter)[0].prefix;
900: newAtlas->updatePoint(*c_iter, &defaultIdx);
901: newSize += defaultIdx.prefix;
902: }
903: section->setAtlas(newAtlas);
904: this->_factory->orderPatch(section, this->getSieve());
905: // Copy over existing values
906: typedef typename alloc_type::template rebind<typename Section_::value_type>::other value_alloc_type;
907: value_alloc_type value_allocator;
908: typename Section_::value_type *newArray = value_allocator.allocate(newSize);
909: for(int i = 0; i < newSize; ++i) {value_allocator.construct(newArray+i, typename Section_::value_type());}
910: ///typename Section_::value_type *newArray = new typename Section_::value_type[newSize];
911: const typename Section_::value_type *array = section->restrict();
913: for(typename Section_::atlas_type::chart_type::const_iterator c_iter = oldChart.begin(); c_iter != oldChart.end(); ++c_iter) {
914: const int& dim = section->getFiberDimension(*c_iter);
915: const int& offset = atlas->restrictPoint(*c_iter)->index;
916: const int& newOffset = newAtlas->restrictPoint(*c_iter)->index;
918: for(int i = 0; i < dim; ++i) {
919: newArray[newOffset+i] = array[offset+i];
920: }
921: }
922: section->replaceStorage(newArray);
923: };
924: public: // Overlap
925: template<typename Sequence>
926: void constructOverlap(const Obj<Sequence>& points, const Obj<send_overlap_type>& sendOverlap, const Obj<recv_overlap_type>& recvOverlap) {
927: point_type *sendBuf = new point_type[points->size()];
928: int size = 0;
929: for(typename Sequence::iterator l_iter = points->begin(); l_iter != points->end(); ++l_iter) {
930: sendBuf[size++] = *l_iter;
931: }
932: int *sizes = new int[this->commSize()]; // The number of points coming from each process
933: int *offsets = new int[this->commSize()+1]; // Prefix sums for sizes
934: int *oldOffs = new int[this->commSize()+1]; // Temporary storage
935: point_type *remotePoints = NULL; // The points from each process
936: int *remoteRanks = NULL; // The rank and number of overlap points of each process that overlaps another
938: // Change to Allgather() for the correct binning algorithm
939: MPI_Gather(&size, 1, MPI_INT, sizes, 1, MPI_INT, 0, this->comm());
940: if (this->commRank() == 0) {
941: offsets[0] = 0;
942: for(int p = 1; p <= this->commSize(); p++) {
943: offsets[p] = offsets[p-1] + sizes[p-1];
944: }
945: remotePoints = new point_type[offsets[this->commSize()]];
946: }
947: MPI_Gatherv(sendBuf, size, MPI_INT, remotePoints, sizes, offsets, MPI_INT, 0, this->comm());
948: delete [] sendBuf;
949: std::map<int, std::map<int, std::set<point_type> > > overlapInfo; // Maps (p,q) to their set of overlap points
951: if (this->commRank() == 0) {
952: for(int p = 0; p < this->commSize(); p++) {
953: std::sort(&remotePoints[offsets[p]], &remotePoints[offsets[p+1]]);
954: }
955: for(int p = 0; p <= this->commSize(); p++) {
956: oldOffs[p] = offsets[p];
957: }
958: for(int p = 0; p < this->commSize(); p++) {
959: for(int q = p+1; q < this->commSize(); q++) {
960: std::set_intersection(&remotePoints[oldOffs[p]], &remotePoints[oldOffs[p+1]],
961: &remotePoints[oldOffs[q]], &remotePoints[oldOffs[q+1]],
962: std::insert_iterator<std::set<point_type> >(overlapInfo[p][q], overlapInfo[p][q].begin()));
963: overlapInfo[q][p] = overlapInfo[p][q];
964: }
965: sizes[p] = overlapInfo[p].size()*2;
966: offsets[p+1] = offsets[p] + sizes[p];
967: }
968: remoteRanks = new int[offsets[this->commSize()]];
969: int k = 0;
970: for(int p = 0; p < this->commSize(); p++) {
971: for(typename std::map<int, std::set<point_type> >::iterator r_iter = overlapInfo[p].begin(); r_iter != overlapInfo[p].end(); ++r_iter) {
972: remoteRanks[k*2] = r_iter->first;
973: remoteRanks[k*2+1] = r_iter->second.size();
974: k++;
975: }
976: }
977: }
978: int numOverlaps; // The number of processes overlapping this process
979: MPI_Scatter(sizes, 1, MPI_INT, &numOverlaps, 1, MPI_INT, 0, this->comm());
980: int *overlapRanks = new int[numOverlaps]; // The rank and overlap size for each overlapping process
981: MPI_Scatterv(remoteRanks, sizes, offsets, MPI_INT, overlapRanks, numOverlaps, MPI_INT, 0, this->comm());
982: point_type *sendPoints = NULL; // The points to send to each process
983: if (this->commRank() == 0) {
984: for(int p = 0, k = 0; p < this->commSize(); p++) {
985: sizes[p] = 0;
986: for(int r = 0; r < (int) overlapInfo[p].size(); r++) {
987: sizes[p] += remoteRanks[k*2+1];
988: k++;
989: }
990: offsets[p+1] = offsets[p] + sizes[p];
991: }
992: sendPoints = new point_type[offsets[this->commSize()]];
993: for(int p = 0, k = 0; p < this->commSize(); p++) {
994: for(typename std::map<int, std::set<point_type> >::iterator r_iter = overlapInfo[p].begin(); r_iter != overlapInfo[p].end(); ++r_iter) {
995: int rank = r_iter->first;
996: for(typename std::set<point_type>::iterator p_iter = (overlapInfo[p][rank]).begin(); p_iter != (overlapInfo[p][rank]).end(); ++p_iter) {
997: sendPoints[k++] = *p_iter;
998: }
999: }
1000: }
1001: }
1002: int numOverlapPoints = 0;
1003: for(int r = 0; r < numOverlaps/2; r++) {
1004: numOverlapPoints += overlapRanks[r*2+1];
1005: }
1006: point_type *overlapPoints = new point_type[numOverlapPoints];
1007: MPI_Scatterv(sendPoints, sizes, offsets, MPI_INT, overlapPoints, numOverlapPoints, MPI_INT, 0, this->comm());
1009: for(int r = 0, k = 0; r < numOverlaps/2; r++) {
1010: int rank = overlapRanks[r*2];
1012: for(int p = 0; p < overlapRanks[r*2+1]; p++) {
1013: point_type point = overlapPoints[k++];
1015: sendOverlap->addArrow(point, rank, point);
1016: recvOverlap->addArrow(rank, point, point);
1017: }
1018: }
1020: delete [] overlapPoints;
1021: delete [] overlapRanks;
1022: delete [] sizes;
1023: delete [] offsets;
1024: delete [] oldOffs;
1025: if (this->commRank() == 0) {
1026: delete [] remoteRanks;
1027: delete [] remotePoints;
1028: delete [] sendPoints;
1029: }
1030: };
1031: void constructOverlap() {
1032: if (this->_calculatedOverlap) return;
1033: this->constructOverlap(this->getSieve()->base(), this->getSendOverlap(), this->getRecvOverlap());
1034: this->constructOverlap(this->getSieve()->cap(), this->getSendOverlap(), this->getRecvOverlap());
1035: if (this->debug()) {
1036: this->_sendOverlap->view("Send overlap");
1037: this->_recvOverlap->view("Receive overlap");
1038: }
1039: this->_calculatedOverlap = true;
1040: };
1041: };
1042: class BoundaryCondition : public ALE::ParallelObject {
1043: public:
1044: typedef double (*function_type)(const double []);
1045: typedef double (*integrator_type)(const double [], const double [], const int, function_type);
1046: protected:
1047: std::string _labelName;
1048: int _marker;
1049: function_type _func;
1050: integrator_type _integrator;
1051: public:
1052: BoundaryCondition(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
1053: ~BoundaryCondition() {};
1054: public:
1055: const std::string& getLabelName() const {return this->_labelName;};
1056: void setLabelName(const std::string& name) {this->_labelName = name;};
1057: int getMarker() const {return this->_marker;};
1058: void setMarker(const int marker) {this->_marker = marker;};
1059: function_type getFunction() const {return this->_func;};
1060: void setFunction(function_type func) {this->_func = func;};
1061: integrator_type getDualIntegrator() const {return this->_integrator;};
1062: void setDualIntegrator(integrator_type integrator) {this->_integrator = integrator;};
1063: public:
1064: double evaluate(const double coords[]) const {return this->_func(coords);};
1065: double integrateDual(const double v0[], const double J[], const int dualIndex) const {return this->_integrator(v0, J, dualIndex, this->_func);};
1066: };
1067: class Discretization : public ALE::ParallelObject {
1068: typedef std::map<std::string, Obj<BoundaryCondition> > boundaryConditions_type;
1069: protected:
1070: boundaryConditions_type _boundaryConditions;
1071: Obj<BoundaryCondition> _exactSolution;
1072: std::map<int,int> _dim2dof;
1073: std::map<int,int> _dim2class;
1074: int _quadSize;
1075: const double *_points;
1076: const double *_weights;
1077: int _basisSize;
1078: const double *_basis;
1079: const double *_basisDer;
1080: const int *_indices;
1081: std::map<int, const int *> _exclusionIndices;
1082: public:
1083: Discretization(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _quadSize(0), _points(NULL), _weights(NULL), _basisSize(0), _basis(NULL), _basisDer(NULL), _indices(NULL) {};
1084: virtual ~Discretization() {
1085: if (this->_indices) {delete [] this->_indices;}
1086: for(std::map<int, const int *>::iterator i_iter = _exclusionIndices.begin(); i_iter != _exclusionIndices.end(); ++i_iter) {
1087: delete [] i_iter->second;
1088: }
1089: };
1090: public:
1091: const bool hasBoundaryCondition() {return (this->_boundaryConditions.find("default") != this->_boundaryConditions.end());};
1092: const Obj<BoundaryCondition>& getBoundaryCondition() {return this->getBoundaryCondition("default");};
1093: void setBoundaryCondition(const Obj<BoundaryCondition>& boundaryCondition) {this->setBoundaryCondition("default", boundaryCondition);};
1094: const Obj<BoundaryCondition>& getBoundaryCondition(const std::string& name) {return this->_boundaryConditions[name];};
1095: void setBoundaryCondition(const std::string& name, const Obj<BoundaryCondition>& boundaryCondition) {this->_boundaryConditions[name] = boundaryCondition;};
1096: Obj<std::set<std::string> > getBoundaryConditions() const {
1097: Obj<std::set<std::string> > names = std::set<std::string>();
1099: for(boundaryConditions_type::const_iterator d_iter = this->_boundaryConditions.begin(); d_iter != this->_boundaryConditions.end(); ++d_iter) {
1100: names->insert(d_iter->first);
1101: }
1102: return names;
1103: };
1104: const Obj<BoundaryCondition>& getExactSolution() {return this->_exactSolution;};
1105: void setExactSolution(const Obj<BoundaryCondition>& exactSolution) {this->_exactSolution = exactSolution;};
1106: const int getQuadratureSize() {return this->_quadSize;};
1107: void setQuadratureSize(const int size) {this->_quadSize = size;};
1108: const double *getQuadraturePoints() {return this->_points;};
1109: void setQuadraturePoints(const double *points) {this->_points = points;};
1110: const double *getQuadratureWeights() {return this->_weights;};
1111: void setQuadratureWeights(const double *weights) {this->_weights = weights;};
1112: const int getBasisSize() {return this->_basisSize;};
1113: void setBasisSize(const int size) {this->_basisSize = size;};
1114: const double *getBasis() {return this->_basis;};
1115: void setBasis(const double *basis) {this->_basis = basis;};
1116: const double *getBasisDerivatives() {return this->_basisDer;};
1117: void setBasisDerivatives(const double *basisDer) {this->_basisDer = basisDer;};
1118: int getNumDof(const int dim) {return this->_dim2dof[dim];};
1119: void setNumDof(const int dim, const int numDof) {this->_dim2dof[dim] = numDof;};
1120: int getDofClass(const int dim) {return this->_dim2class[dim];};
1121: void setDofClass(const int dim, const int dofClass) {this->_dim2class[dim] = dofClass;};
1122: public:
1123: const int *getIndices() {return this->_indices;};
1124: const int *getIndices(const int marker) {
1125: if (!marker) return this->getIndices();
1126: return this->_exclusionIndices[marker];
1127: };
1128: void setIndices(const int *indices) {this->_indices = indices;};
1129: void setIndices(const int *indices, const int marker) {
1130: if (!marker) this->_indices = indices;
1131: this->_exclusionIndices[marker] = indices;
1132: };
1133: template<typename Bundle>
1134: int size(const Obj<Bundle>& mesh) {
1135: const Obj<typename Bundle::label_sequence>& cells = mesh->heightStratum(0);
1136: const Obj<typename Bundle::coneArray> closure = ALE::SieveAlg<Bundle>::closure(mesh, *cells->begin());
1137: const typename Bundle::coneArray::iterator end = closure->end();
1138: int size = 0;
1140: for(typename Bundle::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
1141: size += this->_dim2dof[mesh->depth(*cl_iter)];
1142: }
1143: return size;
1144: };
1145: };
1146: }
1148: namespace ALE {
1149: #ifdef NEW_SECTION
1150: class Mesh : public Bundle<ALE::Sieve<int,int,int>, GeneralSection<int, double> > {
1151: #else
1152: class Mesh : public Bundle<ALE::Sieve<int,int,int> > {
1153: #endif
1154: public:
1155: #ifdef NEW_SECTION
1156: typedef Bundle<ALE::Sieve<int,int,int>, GeneralSection<int, double> > base_type;
1157: #else
1158: typedef Bundle<ALE::Sieve<int,int,int> > base_type;
1159: #endif
1160: typedef base_type::sieve_type sieve_type;
1161: typedef sieve_type::point_type point_type;
1162: typedef base_type::alloc_type alloc_type;
1163: typedef base_type::label_sequence label_sequence;
1164: typedef base_type::real_section_type real_section_type;
1165: typedef base_type::numbering_type numbering_type;
1166: typedef base_type::order_type order_type;
1167: typedef base_type::send_overlap_type send_overlap_type;
1168: typedef base_type::recv_overlap_type recv_overlap_type;
1169: typedef base_type::sieve_alg_type sieve_alg_type;
1170: typedef std::set<std::string> names_type;
1171: typedef std::map<std::string, Obj<Discretization> > discretizations_type;
1172: typedef std::vector<PETSc::Point<3> > holes_type;
1173: protected:
1174: int _dim;
1175: discretizations_type _discretizations;
1176: std::map<int,double> _periodicity;
1177: holes_type _holes;
1178: public:
1179: Mesh(MPI_Comm comm, int dim, int debug = 0) : base_type(comm, debug), _dim(dim) {
1180: ///this->_factory = NumberingFactory::singleton(debug);
1181: };
1182: public: // Accessors
1183: int getDimension() const {return this->_dim;};
1184: void setDimension(const int dim) {this->_dim = dim;};
1185: const Obj<Discretization>& getDiscretization() {return this->getDiscretization("default");};
1186: const Obj<Discretization>& getDiscretization(const std::string& name) {return this->_discretizations[name];};
1187: void setDiscretization(const Obj<Discretization>& disc) {this->setDiscretization("default", disc);};
1188: void setDiscretization(const std::string& name, const Obj<Discretization>& disc) {this->_discretizations[name] = disc;};
1189: Obj<names_type> getDiscretizations() const {
1190: Obj<names_type> names = names_type();
1192: for(discretizations_type::const_iterator d_iter = this->_discretizations.begin(); d_iter != this->_discretizations.end(); ++d_iter) {
1193: names->insert(d_iter->first);
1194: }
1195: return names;
1196: };
1197: bool getPeriodicity(const int d) {return this->_periodicity[d];};
1198: void setPeriodicity(const int d, const double length) {this->_periodicity[d] = length;};
1199: const holes_type& getHoles() const {return this->_holes;};
1200: void addHole(const double hole[]) {
1201: this->_holes.push_back(hole);
1202: };
1203: void copyHoles(const Obj<Mesh>& m) {
1204: const holes_type& holes = m->getHoles();
1206: for(holes_type::const_iterator h_iter = holes.begin(); h_iter != holes.end(); ++h_iter) {
1207: this->_holes.push_back(*h_iter);
1208: }
1209: };
1210: void copy(const Obj<Mesh>& m) {
1211: this->setSieve(m->getSieve());
1212: this->setLabel("height", m->getLabel("height"));
1213: this->_maxHeight = m->height();
1214: this->setLabel("depth", m->getLabel("depth"));
1215: this->_maxDepth = m->depth();
1216: this->setLabel("marker", m->getLabel("marker"));
1217: this->setRealSection("coordinates", m->getRealSection("coordinates"));
1218: this->setArrowSection("orientation", m->getArrowSection("orientation"));
1219: };
1220: public: // Mesh geometry
1221: void computeTriangleGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1222: const double *coords = this->restrict(coordinates, e);
1223: const int dim = 2;
1224: double invDet;
1226: if (v0) {
1227: for(int d = 0; d < dim; d++) {
1228: v0[d] = coords[d];
1229: }
1230: }
1231: if (J) {
1232: for(int d = 0; d < dim; d++) {
1233: for(int f = 0; f < dim; f++) {
1234: J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
1235: }
1236: }
1237: detJ = J[0]*J[3] - J[1]*J[2];
1238: if (detJ < 0.0) {
1239: const double xLength = this->_periodicity[0];
1241: if (xLength != 0.0) {
1242: double v0x = coords[0*dim+0];
1244: if (v0x == 0.0) {
1245: v0x = v0[0] = xLength;
1246: }
1247: for(int f = 0; f < dim; f++) {
1248: const double px = coords[(f+1)*dim+0] == 0.0 ? xLength : coords[(f+1)*dim+0];
1250: J[0*dim+f] = 0.5*(px - v0x);
1251: }
1252: }
1253: detJ = J[0]*J[3] - J[1]*J[2];
1254: }
1255: PetscLogFlopsNoCheck(8 + 3);
1256: }
1257: if (invJ) {
1258: invDet = 1.0/detJ;
1259: invJ[0] = invDet*J[3];
1260: invJ[1] = -invDet*J[1];
1261: invJ[2] = -invDet*J[2];
1262: invJ[3] = invDet*J[0];
1263: PetscLogFlopsNoCheck(5);
1264: }
1265: };
1266: void computeQuadrilateralGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double point[], double v0[], double J[], double invJ[], double& detJ) {
1267: const double *coords = this->restrict(coordinates, e);
1268: const int dim = 2;
1269: double invDet;
1271: if (v0) {
1272: for(int d = 0; d < dim; d++) {
1273: v0[d] = coords[d];
1274: }
1275: }
1276: if (J) {
1277: double x_1 = coords[2] - coords[0];
1278: double y_1 = coords[3] - coords[1];
1279: double x_2 = coords[6] - coords[0];
1280: double y_2 = coords[7] - coords[1];
1281: double x_3 = coords[4] - coords[0];
1282: double y_3 = coords[5] - coords[1];
1284: J[0] = x_1 + (x_3 - x_1 - x_2)*point[1];
1285: J[1] = x_2 + (x_3 - x_1 - x_2)*point[0];
1286: J[2] = y_1 + (y_3 - y_1 - y_2)*point[1];
1287: J[3] = y_1 + (y_3 - y_1 - y_2)*point[0];
1288: detJ = J[0]*J[3] - J[1]*J[2];
1289: PetscLogFlopsNoCheck(6 + 16 + 3);
1290: }
1291: if (invJ) {
1292: invDet = 1.0/detJ;
1293: invJ[0] = invDet*J[3];
1294: invJ[1] = -invDet*J[1];
1295: invJ[2] = -invDet*J[2];
1296: invJ[3] = invDet*J[0];
1297: PetscLogFlopsNoCheck(5);
1298: }
1299: };
1300: void computeTetrahedronGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1301: const double *coords = this->restrict(coordinates, e);
1302: const int dim = 3;
1303: double invDet;
1305: if (v0) {
1306: for(int d = 0; d < dim; d++) {
1307: v0[d] = coords[d];
1308: }
1309: }
1310: if (J) {
1311: for(int d = 0; d < dim; d++) {
1312: for(int f = 0; f < dim; f++) {
1313: J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
1314: }
1315: }
1316: // The minus sign is here since I orient the first face to get the outward normal
1317: detJ = -(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
1318: J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
1319: J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
1320: PetscLogFlopsNoCheck(18 + 12);
1321: }
1322: if (invJ) {
1323: invDet = -1.0/detJ;
1324: invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
1325: invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
1326: invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
1327: invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
1328: invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
1329: invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
1330: invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
1331: invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
1332: invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
1333: PetscLogFlopsNoCheck(37);
1334: }
1335: };
1336: void computeHexahedralGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double point[], double v0[], double J[], double invJ[], double& detJ) {
1337: const double *coords = this->restrict(coordinates, e);
1338: const int dim = 3;
1339: double invDet;
1341: if (v0) {
1342: for(int d = 0; d < dim; d++) {
1343: v0[d] = coords[d];
1344: }
1345: }
1346: if (J) {
1347: double x_1 = coords[3] - coords[0];
1348: double y_1 = coords[4] - coords[1];
1349: double z_1 = coords[5] - coords[2];
1350: double x_2 = coords[6] - coords[0];
1351: double y_2 = coords[7] - coords[1];
1352: double z_2 = coords[8] - coords[2];
1353: double x_3 = coords[9] - coords[0];
1354: double y_3 = coords[10] - coords[1];
1355: double z_3 = coords[11] - coords[2];
1356: double x_4 = coords[12] - coords[0];
1357: double y_4 = coords[13] - coords[1];
1358: double z_4 = coords[14] - coords[2];
1359: double x_5 = coords[15] - coords[0];
1360: double y_5 = coords[16] - coords[1];
1361: double z_5 = coords[17] - coords[2];
1362: double x_6 = coords[18] - coords[0];
1363: double y_6 = coords[19] - coords[1];
1364: double z_6 = coords[20] - coords[2];
1365: double x_7 = coords[21] - coords[0];
1366: double y_7 = coords[22] - coords[1];
1367: double z_7 = coords[23] - coords[2];
1368: double g_x = x_1 - x_2 + x_3 + x_4 - x_5 + x_6 - x_7;
1369: double g_y = y_1 - y_2 + y_3 + y_4 - y_5 + y_6 - y_7;
1370: double g_z = z_1 - z_2 + z_3 + z_4 - z_5 + z_6 - z_7;
1372: 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];
1373: 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];
1374: 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];
1375: 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];
1376: 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];
1377: 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];
1378: 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];
1379: 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];
1380: 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];
1381: detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
1382: J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
1383: J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
1384: PetscLogFlopsNoCheck(39 + 81 + 12);
1385: }
1386: if (invJ) {
1387: invDet = 1.0/detJ;
1388: invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
1389: invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
1390: invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
1391: invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
1392: invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
1393: invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
1394: invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
1395: invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
1396: invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
1397: PetscLogFlopsNoCheck(37);
1398: }
1399: };
1400: void computeElementGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
1401: if (this->_dim == 2) {
1402: computeTriangleGeometry(coordinates, e, v0, J, invJ, detJ);
1403: } else if (this->_dim == 3) {
1404: computeTetrahedronGeometry(coordinates, e, v0, J, invJ, detJ);
1405: } else {
1406: throw ALE::Exception("Unsupport dimension for element geometry computation");
1407: }
1408: };
1409: void computeLineFaceGeometry(const point_type& cell, const point_type& face, const int f, const double cellInvJ[], double invJ[], double& detJ, double normal[], double tangent[]) {
1410: const arrow_section_type::point_type arrow(cell, face);
1411: const bool reversed = (this->getArrowSection("orientation")->restrictPoint(arrow)[0] == -2);
1412: const int dim = this->getDimension();
1413: double norm = 0.0;
1414: double *vec = tangent;
1416: if (f == 0) {
1417: vec[0] = 0.0; vec[1] = -1.0;
1418: } else if (f == 1) {
1419: vec[0] = 0.70710678; vec[1] = 0.70710678;
1420: } else if (f == 2) {
1421: vec[0] = -1.0; vec[1] = 0.0;
1422: }
1423: for(int d = 0; d < dim; ++d) {
1424: normal[d] = 0.0;
1425: for(int e = 0; e < dim; ++e) normal[d] += cellInvJ[e*dim+d]*vec[e];
1426: if (reversed) normal[d] = -normal[d];
1427: norm += normal[d]*normal[d];
1428: }
1429: norm = std::sqrt(norm);
1430: for(int d = 0; d < dim; ++d) {
1431: normal[d] /= norm;
1432: }
1433: tangent[0] = normal[1];
1434: tangent[1] = -normal[0];
1435: if (this->debug()) {
1436: std::cout << "Cell: " << cell << " Face: " << face << "("<<f<<")" << std::endl;
1437: for(int d = 0; d < dim; ++d) {
1438: std::cout << "Normal["<<d<<"]: " << normal[d] << " Tangent["<<d<<"]: " << tangent[d] << std::endl;
1439: }
1440: }
1441: // Now get 1D Jacobian info
1442: // Should be a way to get this directly
1443: const double *coords = this->restrict(this->getRealSection("coordinates"), face);
1444: detJ = std::sqrt(PetscSqr(coords[1*2+0] - coords[0*2+0]) + PetscSqr(coords[1*2+1] - coords[0*2+1]))/2.0;
1445: invJ[0] = 1.0/detJ;
1446: };
1447: void computeTriangleFaceGeometry(const point_type& cell, const point_type& face, const int f, const double cellInvJ[], double invJ[], double& detJ, double normal[], double tangent[]) {
1448: const arrow_section_type::point_type arrow(cell, face);
1449: const bool reversed = this->getArrowSection("orientation")->restrictPoint(arrow)[0] < 0;
1450: const int dim = this->getDimension();
1451: const int faceDim = dim-1;
1452: double norm = 0.0;
1453: double *vec = tangent;
1455: if (f == 0) {
1456: vec[0] = 0.0; vec[1] = 0.0; vec[2] = -1.0;
1457: } else if (f == 1) {
1458: vec[0] = 0.0; vec[1] = -1.0; vec[2] = 0.0;
1459: } else if (f == 2) {
1460: vec[0] = 0.57735027; vec[1] = 0.57735027; vec[2] = 0.57735027;
1461: } else if (f == 3) {
1462: vec[0] = -1.0; vec[1] = 0.0; vec[2] = 0.0;
1463: }
1464: for(int d = 0; d < dim; ++d) {
1465: normal[d] = 0.0;
1466: for(int e = 0; e < dim; ++e) normal[d] += cellInvJ[e*dim+d]*vec[e];
1467: if (reversed) normal[d] = -normal[d];
1468: norm += normal[d]*normal[d];
1469: }
1470: norm = std::sqrt(norm);
1471: for(int d = 0; d < dim; ++d) {
1472: normal[d] /= norm;
1473: }
1474: // Get tangents
1475: tangent[0] = normal[1] - normal[2];
1476: tangent[1] = normal[2] - normal[0];
1477: tangent[2] = normal[0] - normal[1];
1478: norm = 0.0;
1479: for(int d = 0; d < dim; ++d) {
1480: norm += tangent[d]*tangent[d];
1481: }
1482: norm = std::sqrt(norm);
1483: for(int d = 0; d < dim; ++d) {
1484: tangent[d] /= norm;
1485: }
1486: tangent[3] = normal[1]*tangent[2] - normal[2]*tangent[1];
1487: tangent[4] = normal[2]*tangent[0] - normal[0]*tangent[2];
1488: tangent[5] = normal[0]*tangent[1] - normal[1]*tangent[0];
1489: if (this->debug()) {
1490: std::cout << "Cell: " << cell << " Face: " << face << "("<<f<<")" << std::endl;
1491: for(int d = 0; d < dim; ++d) {
1492: std::cout << "Normal["<<d<<"]: " << normal[d] << " TangentA["<<d<<"]: " << tangent[d] << " TangentB["<<d<<"]: " << tangent[dim+d] << std::endl;
1493: }
1494: }
1495: // Now get 2D Jacobian info
1496: // Should be a way to get this directly
1497: const double *coords = this->restrict(this->getRealSection("coordinates"), face);
1498: // Rotate so that normal in z
1499: double invR[9], R[9];
1500: double detR, invDetR;
1501: for(int d = 0; d < dim; d++) {
1502: invR[d*dim+0] = tangent[d];
1503: invR[d*dim+1] = tangent[dim+d];
1504: invR[d*dim+2] = normal[d];
1505: }
1506: invDetR = (invR[0*3+0]*(invR[1*3+1]*invR[2*3+2] - invR[1*3+2]*invR[2*3+1]) +
1507: invR[0*3+1]*(invR[1*3+2]*invR[2*3+0] - invR[1*3+0]*invR[2*3+2]) +
1508: invR[0*3+2]*(invR[1*3+0]*invR[2*3+1] - invR[1*3+1]*invR[2*3+0]));
1509: detR = 1.0/invDetR;
1510: R[0*3+0] = detR*(invR[1*3+1]*invR[2*3+2] - invR[1*3+2]*invR[2*3+1]);
1511: R[0*3+1] = detR*(invR[0*3+2]*invR[2*3+1] - invR[0*3+1]*invR[2*3+2]);
1512: R[0*3+2] = detR*(invR[0*3+1]*invR[1*3+2] - invR[0*3+2]*invR[1*3+1]);
1513: R[1*3+0] = detR*(invR[1*3+2]*invR[2*3+0] - invR[1*3+0]*invR[2*3+2]);
1514: R[1*3+1] = detR*(invR[0*3+0]*invR[2*3+2] - invR[0*3+2]*invR[2*3+0]);
1515: R[1*3+2] = detR*(invR[0*3+2]*invR[1*3+0] - invR[0*3+0]*invR[1*3+2]);
1516: R[2*3+0] = detR*(invR[1*3+0]*invR[2*3+1] - invR[1*3+1]*invR[2*3+0]);
1517: R[2*3+1] = detR*(invR[0*3+1]*invR[2*3+0] - invR[0*3+0]*invR[2*3+1]);
1518: R[2*3+2] = detR*(invR[0*3+0]*invR[1*3+1] - invR[0*3+1]*invR[1*3+0]);
1519: for(int d = 0; d < dim; d++) {
1520: for(int e = 0; e < dim; e++) {
1521: invR[d*dim+e] = 0.0;
1522: for(int g = 0; g < dim; g++) {
1523: invR[d*dim+e] += R[e*dim+g]*coords[d*dim+g];
1524: }
1525: }
1526: }
1527: for(int d = dim-1; d >= 0; --d) {
1528: invR[d*dim+2] -= invR[0*dim+2];
1529: if (this->debug() && (d == dim-1)) {
1530: double ref[9];
1531: for(int q = 0; q < dim; q++) {
1532: for(int e = 0; e < dim; e++) {
1533: ref[q*dim+e] = 0.0;
1534: for(int g = 0; g < dim; g++) {
1535: ref[q*dim+e] += cellInvJ[e*dim+g]*coords[q*dim+g];
1536: }
1537: }
1538: }
1539: std::cout << "f: " << f << std::endl;
1540: std::cout << this->printMatrix(std::string("coords"), dim, dim, coords) << std::endl;
1541: std::cout << this->printMatrix(std::string("ref coords"), dim, dim, ref) << std::endl;
1542: std::cout << this->printMatrix(std::string("R"), dim, dim, R) << std::endl;
1543: std::cout << this->printMatrix(std::string("invR"), dim, dim, invR) << std::endl;
1544: }
1545: if (fabs(invR[d*dim+2]) > 1.0e-8) {
1546: throw ALE::Exception("Invalid rotation");
1547: }
1548: }
1549: double J[4];
1550: for(int d = 0; d < faceDim; d++) {
1551: for(int e = 0; e < faceDim; e++) {
1552: J[d*faceDim+e] = 0.5*(invR[(e+1)*dim+d] - invR[0*dim+d]);
1553: }
1554: }
1555: detJ = fabs(J[0]*J[3] - J[1]*J[2]);
1556: // Probably need something here if detJ < 0
1557: const double invDet = 1.0/detJ;
1558: invJ[0] = invDet*J[3];
1559: invJ[1] = -invDet*J[1];
1560: invJ[2] = -invDet*J[2];
1561: invJ[3] = invDet*J[0];
1562: };
1563: void computeFaceGeometry(const point_type& cell, const point_type& face, const int f, const double cellInvJ[], double invJ[], double& detJ, double normal[], double tangent[]) {
1564: if (this->_dim == 2) {
1565: computeLineFaceGeometry(cell, face, f, cellInvJ, invJ, detJ, normal, tangent);
1566: } else if (this->_dim == 3) {
1567: computeTriangleFaceGeometry(cell, face, f, cellInvJ, invJ, detJ, normal, tangent);
1568: } else {
1569: throw ALE::Exception("Unsupport dimension for element geometry computation");
1570: }
1571: };
1572: double getMaxVolume() {
1573: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1574: const Obj<label_sequence>& cells = this->heightStratum(0);
1575: const int dim = this->getDimension();
1576: double v0[3], J[9], invJ[9], detJ, refVolume = 0.0, maxVolume = 0.0;
1578: if (dim == 1) refVolume = 2.0;
1579: if (dim == 2) refVolume = 2.0;
1580: if (dim == 3) refVolume = 4.0/3.0;
1581: for(label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
1582: this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
1583: maxVolume = std::max(maxVolume, detJ*refVolume);
1584: }
1585: return maxVolume;
1586: };
1587: // Find the cell in which this point lies (stupid algorithm)
1588: point_type locatePoint_2D(const real_section_type::value_type point[]) {
1589: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1590: const Obj<label_sequence>& cells = this->heightStratum(0);
1591: const int embedDim = 2;
1592: double v0[2], J[4], invJ[4], detJ;
1594: for(label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
1595: this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
1596: double xi = invJ[0*embedDim+0]*(point[0] - v0[0]) + invJ[0*embedDim+1]*(point[1] - v0[1]);
1597: double eta = invJ[1*embedDim+0]*(point[0] - v0[0]) + invJ[1*embedDim+1]*(point[1] - v0[1]);
1599: if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) {
1600: return *c_iter;
1601: }
1602: }
1603: throw ALE::Exception("Could not locate point");
1604: };
1605: // Assume a simplex and 3D
1606: point_type locatePoint_3D(const real_section_type::value_type point[]) {
1607: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1608: const Obj<label_sequence>& cells = this->heightStratum(0);
1609: const int embedDim = 3;
1610: double v0[3], J[9], invJ[9], detJ;
1612: for(label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
1613: this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
1614: 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]);
1615: 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]);
1616: 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]);
1618: if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) {
1619: return *c_iter;
1620: }
1621: }
1622: throw ALE::Exception("Could not locate point");
1623: };
1624: point_type locatePoint(const real_section_type::value_type point[], point_type guess = -1) {
1625: //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
1626: if (guess != -1) {
1627: return guess;
1628: }else if (this->_dim == 2) {
1629: return locatePoint_2D(point);
1630: } else if (this->_dim == 3) {
1631: return locatePoint_3D(point);
1632: } else {
1633: throw ALE::Exception("No point location for mesh dimension");
1634: }
1635: };
1636: public: // Discretization
1637: void markBoundaryCells(const std::string& name, const int marker = 1, const int newMarker = 2, const bool onlyVertices = false) {
1638: const Obj<label_type>& label = this->getLabel(name);
1639: const Obj<label_sequence>& boundary = this->getLabelStratum(name, marker);
1640: const Obj<sieve_type>& sieve = this->getSieve();
1642: if (!onlyVertices) {
1643: const label_sequence::iterator end = boundary->end();
1645: for(label_sequence::iterator e_iter = boundary->begin(); e_iter != end; ++e_iter) {
1646: if (this->height(*e_iter) == 1) {
1647: const point_type cell = *sieve->support(*e_iter)->begin();
1649: this->setValue(label, cell, newMarker);
1650: }
1651: }
1652: } else {
1653: const label_sequence::iterator end = boundary->end();
1654: const int depth = this->depth();
1656: for(label_sequence::iterator v_iter = boundary->begin(); v_iter != end; ++v_iter) {
1657: const Obj<supportArray> support = sieve->nSupport(*v_iter, depth);
1658: const supportArray::iterator sEnd = support->end();
1660: for(supportArray::iterator c_iter = support->begin(); c_iter != sEnd; ++c_iter) {
1661: const Obj<sieve_type::traits::coneSequence>& cone = sieve->cone(*c_iter);
1662: const sieve_type::traits::coneSequence::iterator cEnd = cone->end();
1664: for(sieve_type::traits::coneSequence::iterator e_iter = cone->begin(); e_iter != cEnd; ++e_iter) {
1665: if (sieve->support(*e_iter)->size() == 1) {
1666: this->setValue(label, *c_iter, newMarker);
1667: break;
1668: }
1669: }
1670: }
1671: }
1672: }
1673: };
1674: int setFiberDimensions(const Obj<real_section_type>& s, const Obj<names_type>& discs, names_type& bcLabels) {
1675: const int debug = this->debug();
1676: int maxDof = 0;
1678: for(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter) {
1679: s->addSpace();
1680: }
1681: for(int d = 0; d <= this->_dim; ++d) {
1682: int numDof = 0;
1683: int f = 0;
1685: for(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
1686: const Obj<ALE::Discretization>& disc = this->getDiscretization(*f_iter);
1688: numDof += disc->getNumDof(d);
1689: s->setFiberDimension(this->depthStratum(d), disc->getNumDof(d), f);
1690: }
1691: s->setFiberDimension(this->depthStratum(d), numDof);
1692: maxDof = std::max(maxDof, numDof);
1693: }
1694: // Process exclusions
1695: int f = 0;
1697: for(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
1698: const Obj<ALE::Discretization>& disc = this->getDiscretization(*f_iter);
1699: std::string labelName = "exclude-"+*f_iter;
1700: std::set<point_type> seen;
1702: if (this->hasLabel(labelName)) {
1703: const Obj<label_type>& label = this->getLabel(labelName);
1704: const Obj<label_sequence>& exclusion = this->getLabelStratum(labelName, 1);
1705: const label_sequence::iterator end = exclusion->end();
1706: if (debug > 1) {label->view(labelName.c_str());}
1708: for(label_sequence::iterator e_iter = exclusion->begin(); e_iter != end; ++e_iter) {
1709: const Obj<coneArray> closure = ALE::SieveAlg<ALE::Mesh>::closure(this, this->getArrowSection("orientation"), *e_iter);
1710: const coneArray::iterator cEnd = closure->end();
1712: for(coneArray::iterator c_iter = closure->begin(); c_iter != cEnd; ++c_iter) {
1713: if (seen.find(*c_iter) != seen.end()) continue;
1714: if (this->getValue(label, *c_iter) == 1) {
1715: seen.insert(*c_iter);
1716: s->setFiberDimension(*c_iter, 0, f);
1717: s->addFiberDimension(*c_iter, -disc->getNumDof(this->depth(*c_iter)));
1718: if (debug > 1) {std::cout << " cell: " << *c_iter << " dim: " << disc->getNumDof(this->depth(*c_iter)) << std::endl;}
1719: }
1720: }
1721: }
1722: }
1723: }
1724: // Process constraints
1725: f = 0;
1726: for(std::set<std::string>::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
1727: const Obj<ALE::Discretization>& disc = this->getDiscretization(*f_iter);
1728: const Obj<std::set<std::string> > bcs = disc->getBoundaryConditions();
1729: std::string excludeName = "exclude-"+*f_iter;
1731: for(std::set<std::string>::const_iterator bc_iter = bcs->begin(); bc_iter != bcs->end(); ++bc_iter) {
1732: const Obj<ALE::BoundaryCondition>& bc = disc->getBoundaryCondition(*bc_iter);
1733: const Obj<label_sequence>& boundary = this->getLabelStratum(bc->getLabelName(), bc->getMarker());
1735: bcLabels.insert(bc->getLabelName());
1736: if (this->hasLabel(excludeName)) {
1737: const Obj<label_type>& label = this->getLabel(excludeName);
1739: for(label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
1740: if (!this->getValue(label, *e_iter)) {
1741: s->addConstraintDimension(*e_iter, disc->getNumDof(this->depth(*e_iter)));
1742: s->setConstraintDimension(*e_iter, disc->getNumDof(this->depth(*e_iter)), f);
1743: }
1744: }
1745: } else {
1746: for(label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
1747: s->addConstraintDimension(*e_iter, disc->getNumDof(this->depth(*e_iter)));
1748: s->setConstraintDimension(*e_iter, disc->getNumDof(this->depth(*e_iter)), f);
1749: }
1750: }
1751: }
1752: }
1753: return maxDof;
1754: };
1755: void calculateIndices() {
1756: // Should have an iterator over the whole tree
1757: Obj<names_type> discs = this->getDiscretizations();
1758: Obj<Mesh> mesh = this;
1759: const int debug = this->debug();
1760: std::map<std::string, std::pair<int, int*> > indices;
1762: mesh.addRef();
1763: for(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1764: const Obj<Discretization>& disc = this->getDiscretization(*d_iter);
1766: indices[*d_iter] = std::pair<int, int*>(0, new int[disc->size(mesh)]);
1767: disc->setIndices(indices[*d_iter].second);
1768: }
1769: const Obj<label_sequence>& cells = this->heightStratum(0);
1770: const Obj<coneArray> closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), *cells->begin());
1771: const coneArray::iterator end = closure->end();
1772: int offset = 0;
1774: if (debug > 1) {std::cout << "Closure for first element" << std::endl;}
1775: for(coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
1776: const int dim = this->depth(*cl_iter);
1778: if (debug > 1) {std::cout << " point " << *cl_iter << " depth " << dim << std::endl;}
1779: for(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1780: const Obj<Discretization>& disc = this->getDiscretization(*d_iter);
1781: const int num = disc->getNumDof(dim);
1783: if (debug > 1) {std::cout << " disc " << disc->getName() << " numDof " << num << std::endl;}
1784: for(int o = 0; o < num; ++o) {
1785: indices[*d_iter].second[indices[*d_iter].first++] = offset++;
1786: }
1787: }
1788: }
1789: if (debug > 1) {
1790: for(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1791: const Obj<Discretization>& disc = this->getDiscretization(*d_iter);
1793: std::cout << "Discretization " << disc->getName() << " indices:";
1794: for(int i = 0; i < indices[*d_iter].first; ++i) {
1795: std::cout << " " << indices[*d_iter].second[i];
1796: }
1797: std::cout << std::endl;
1798: }
1799: }
1800: };
1801: void calculateIndicesExcluded(const Obj<real_section_type>& s, const Obj<names_type>& discs) {
1802: typedef std::map<std::string, std::pair<int, indexSet> > indices_type;
1803: const Obj<label_type>& indexLabel = this->createLabel("cellExclusion");
1804: Obj<Mesh> mesh = this;
1805: const int debug = this->debug();
1806: int marker = 0;
1807: std::map<indices_type, int> indexMap;
1808: indices_type indices;
1810: mesh.addRef();
1811: for(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
1812: const Obj<Discretization>& disc = this->getDiscretization(*d_iter);
1813: const int size = disc->size(mesh);
1815: indices[*d_iter].second.resize(size);
1816: }
1817: const names_type::const_iterator dBegin = discs->begin();
1818: const names_type::const_iterator dEnd = discs->end();
1819: std::set<point_type> seen;
1820: int f = 0;
1822: for(names_type::const_iterator f_iter = dBegin; f_iter != dEnd; ++f_iter, ++f) {
1823: std::string labelName = "exclude-"+*f_iter;
1825: if (this->hasLabel(labelName)) {
1826: const Obj<label_sequence>& exclusion = this->getLabelStratum(labelName, 1);
1827: const label_sequence::iterator end = exclusion->end();
1829: if (debug > 1) {std::cout << "Processing exclusion " << labelName << std::endl;}
1830: for(label_sequence::iterator e_iter = exclusion->begin(); e_iter != end; ++e_iter) {
1831: if (this->height(*e_iter)) continue;
1832: const Obj<coneArray> closure = ALE::SieveAlg<ALE::Mesh>::closure(this, this->getArrowSection("orientation"), *e_iter);
1833: const coneArray::iterator clEnd = closure->end();
1834: int offset = 0;
1836: if (debug > 1) {std::cout << " Closure for cell " << *e_iter << std::endl;}
1837: for(coneArray::iterator cl_iter = closure->begin(); cl_iter != clEnd; ++cl_iter) {
1838: int g = 0;
1840: if (debug > 1) {std::cout << " point " << *cl_iter << std::endl;}
1841: for(names_type::const_iterator g_iter = dBegin; g_iter != dEnd; ++g_iter, ++g) {
1842: const int fDim = s->getFiberDimension(*cl_iter, g);
1844: if (debug > 1) {std::cout << " disc " << *g_iter << " numDof " << fDim << std::endl;}
1845: for(int d = 0; d < fDim; ++d) {
1846: indices[*g_iter].second[indices[*g_iter].first++] = offset++;
1847: }
1848: }
1849: }
1850: const std::map<indices_type, int>::iterator entry = indexMap.find(indices);
1852: if (debug > 1) {
1853: for(std::map<indices_type, int>::iterator i_iter = indexMap.begin(); i_iter != indexMap.end(); ++i_iter) {
1854: for(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
1855: std::cout << "Discretization (" << i_iter->second << ") " << *g_iter << " indices:";
1856: for(int i = 0; i < ((indices_type) i_iter->first)[*g_iter].first; ++i) {
1857: std::cout << " " << ((indices_type) i_iter->first)[*g_iter].second[i];
1858: }
1859: std::cout << std::endl;
1860: }
1861: std::cout << "Comparison: " << (indices == i_iter->first) << std::endl;
1862: }
1863: for(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
1864: std::cout << "Discretization " << *g_iter << " indices:";
1865: for(int i = 0; i < indices[*g_iter].first; ++i) {
1866: std::cout << " " << indices[*g_iter].second[i];
1867: }
1868: std::cout << std::endl;
1869: }
1870: }
1871: if (entry != indexMap.end()) {
1872: this->setValue(indexLabel, *e_iter, entry->second);
1873: if (debug > 1) {std::cout << " Found existing indices with marker " << entry->second << std::endl;}
1874: } else {
1875: indexMap[indices] = ++marker;
1876: this->setValue(indexLabel, *e_iter, marker);
1877: if (debug > 1) {std::cout << " Created new indices with marker " << marker << std::endl;}
1878: }
1879: for(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
1880: indices[*g_iter].first = 0;
1881: for(unsigned int i = 0; i < indices[*g_iter].second.size(); ++i) indices[*g_iter].second[i] = 0;
1882: }
1883: }
1884: }
1885: }
1886: if (debug > 1) {indexLabel->view("cellExclusion");}
1887: for(std::map<indices_type, int>::iterator i_iter = indexMap.begin(); i_iter != indexMap.end(); ++i_iter) {
1888: if (debug > 1) {std::cout << "Setting indices for marker " << i_iter->second << std::endl;}
1889: for(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
1890: const Obj<Discretization>& disc = this->getDiscretization(*g_iter);
1891: const indexSet indSet = ((indices_type) i_iter->first)[*g_iter].second;
1892: const int size = indSet.size();
1893: int *_indices = new int[size];
1895: if (debug > 1) {std::cout << " field " << *g_iter << std::endl;}
1896: for(int i = 0; i < size; ++i) {
1897: _indices[i] = indSet[i];
1898: if (debug > 1) {std::cout << " indices["<<i<<"] = " << _indices[i] << std::endl;}
1899: }
1900: disc->setIndices(_indices, i_iter->second);
1901: }
1902: }
1903: };
1904: void setupField(const Obj<real_section_type>& s, const int cellMarker = 2, const bool noUpdate = false) {
1905: const Obj<names_type>& discs = this->getDiscretizations();
1906: const int debug = s->debug();
1907: names_type bcLabels;
1908: int maxDof;
1910: maxDof = this->setFiberDimensions(s, discs, bcLabels);
1911: this->calculateIndices();
1912: this->calculateIndicesExcluded(s, discs);
1913: this->allocate(s);
1914: s->defaultConstraintDof();
1915: const Obj<label_type>& cellExclusion = this->getLabel("cellExclusion");
1917: if (debug > 1) {std::cout << "Setting boundary values" << std::endl;}
1918: for(names_type::const_iterator n_iter = bcLabels.begin(); n_iter != bcLabels.end(); ++n_iter) {
1919: const Obj<label_sequence>& boundaryCells = this->getLabelStratum(*n_iter, cellMarker);
1920: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
1921: const Obj<names_type>& discs = this->getDiscretizations();
1922: const point_type firstCell = *boundaryCells->begin();
1923: const int numFields = discs->size();
1924: real_section_type::value_type *values = new real_section_type::value_type[this->sizeWithBC(s, firstCell)];
1925: int *dofs = new int[maxDof];
1926: int *v = new int[numFields];
1927: double *v0 = new double[this->getDimension()];
1928: double *J = new double[this->getDimension()*this->getDimension()];
1929: double detJ;
1931: for(label_sequence::iterator c_iter = boundaryCells->begin(); c_iter != boundaryCells->end(); ++c_iter) {
1932: const Obj<coneArray> closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), *c_iter);
1933: const coneArray::iterator end = closure->end();
1935: if (debug > 1) {std::cout << " Boundary cell " << *c_iter << std::endl;}
1936: this->computeElementGeometry(coordinates, *c_iter, v0, J, PETSC_NULL, detJ);
1937: for(int f = 0; f < numFields; ++f) v[f] = 0;
1938: for(coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
1939: const int cDim = s->getConstraintDimension(*cl_iter);
1940: int off = 0;
1941: int f = 0;
1942: int i = -1;
1944: if (debug > 1) {std::cout << " point " << *cl_iter << std::endl;}
1945: if (cDim) {
1946: if (debug > 1) {std::cout << " constrained excMarker: " << this->getValue(cellExclusion, *c_iter) << std::endl;}
1947: for(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
1948: const Obj<ALE::Discretization>& disc = this->getDiscretization(*f_iter);
1949: const Obj<names_type> bcs = disc->getBoundaryConditions();
1950: const int fDim = s->getFiberDimension(*cl_iter, f);//disc->getNumDof(this->depth(*cl_iter));
1951: const int *indices = disc->getIndices(this->getValue(cellExclusion, *c_iter));
1952: int b = 0;
1954: for(names_type::const_iterator bc_iter = bcs->begin(); bc_iter != bcs->end(); ++bc_iter, ++b) {
1955: const Obj<ALE::BoundaryCondition>& bc = disc->getBoundaryCondition(*bc_iter);
1956: const int value = this->getValue(this->getLabel(bc->getLabelName()), *cl_iter);
1958: if (b > 0) v[f] -= fDim;
1959: if (value == bc->getMarker()) {
1960: if (debug > 1) {std::cout << " field " << *f_iter << " marker " << value << std::endl;}
1961: for(int d = 0; d < fDim; ++d, ++v[f]) {
1962: dofs[++i] = off+d;
1963: if (!noUpdate) values[indices[v[f]]] = (*bc->getDualIntegrator())(v0, J, v[f], bc->getFunction());
1964: if (debug > 1) {std::cout << " setting values["<<indices[v[f]]<<"] = " << values[indices[v[f]]] << std::endl;}
1965: }
1966: // Allow only one condition per point
1967: ++b;
1968: break;
1969: } else {
1970: if (debug > 1) {std::cout << " field " << *f_iter << std::endl;}
1971: for(int d = 0; d < fDim; ++d, ++v[f]) {
1972: values[indices[v[f]]] = 0.0;
1973: if (debug > 1) {std::cout << " setting values["<<indices[v[f]]<<"] = " << values[indices[v[f]]] << std::endl;}
1974: }
1975: }
1976: }
1977: if (b == 0) {
1978: if (debug > 1) {std::cout << " field " << *f_iter << std::endl;}
1979: for(int d = 0; d < fDim; ++d, ++v[f]) {
1980: values[indices[v[f]]] = 0.0;
1981: if (debug > 1) {std::cout << " setting values["<<indices[v[f]]<<"] = " << values[indices[v[f]]] << std::endl;}
1982: }
1983: }
1984: off += fDim;
1985: }
1986: if (i != cDim-1) {throw ALE::Exception("Invalid constraint initialization");}
1987: s->setConstraintDof(*cl_iter, dofs);
1988: } else {
1989: if (debug > 1) {std::cout << " unconstrained" << std::endl;}
1990: for(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
1991: const Obj<ALE::Discretization>& disc = this->getDiscretization(*f_iter);
1992: const int fDim = s->getFiberDimension(*cl_iter, f);//disc->getNumDof(this->depth(*cl_iter));
1993: const int *indices = disc->getIndices(this->getValue(cellExclusion, *c_iter));
1995: if (debug > 1) {std::cout << " field " << *f_iter << std::endl;}
1996: for(int d = 0; d < fDim; ++d, ++v[f]) {
1997: values[indices[v[f]]] = 0.0;
1998: if (debug > 1) {std::cout << " setting values["<<indices[v[f]]<<"] = " << values[indices[v[f]]] << std::endl;}
1999: }
2000: }
2001: }
2002: }
2003: if (debug > 1) {
2004: const Obj<coneArray> closure = sieve_alg_type::closure(this, this->getArrowSection("orientation"), *c_iter);
2005: const coneArray::iterator end = closure->end();
2007: for(int f = 0; f < numFields; ++f) v[f] = 0;
2008: for(coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
2009: int f = 0;
2010: for(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
2011: const Obj<ALE::Discretization>& disc = this->getDiscretization(*f_iter);
2012: const int fDim = s->getFiberDimension(*cl_iter, f);
2013: const int *indices = disc->getIndices(this->getValue(cellExclusion, *c_iter));
2015: for(int d = 0; d < fDim; ++d, ++v[f]) {
2016: std::cout << " "<<*f_iter<<"-value["<<indices[v[f]]<<"] " << values[indices[v[f]]] << std::endl;
2017: }
2018: }
2019: }
2020: }
2021: if (!noUpdate) {
2022: this->updateAll(s, *c_iter, values);
2023: }
2024: }
2025: delete [] dofs;
2026: delete [] values;
2027: delete [] v0;
2028: delete [] J;
2029: }
2030: if (debug > 1) {s->view("");}
2031: };
2032: public:
2033: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
2034: if (comm == MPI_COMM_NULL) {
2035: comm = this->comm();
2036: }
2037: if (name == "") {
2038: PetscPrintf(comm, "viewing a Mesh\n");
2039: } else {
2040: PetscPrintf(comm, "viewing Mesh '%s'\n", name.c_str());
2041: }
2042: this->getSieve()->view("mesh sieve", comm);
2043: Obj<names_type> sections = this->getRealSections();
2045: for(names_type::iterator name = sections->begin(); name != sections->end(); ++name) {
2046: this->getRealSection(*name)->view(*name);
2047: }
2048: sections = this->getIntSections();
2049: for(names_type::iterator name = sections->begin(); name != sections->end(); ++name) {
2050: this->getIntSection(*name)->view(*name);
2051: }
2052: sections = this->getArrowSections();
2053: for(names_type::iterator name = sections->begin(); name != sections->end(); ++name) {
2054: this->getArrowSection(*name)->view(*name);
2055: }
2056: };
2057: template<typename value_type>
2058: static std::string printMatrix(const std::string& name, const int rows, const int cols, const value_type matrix[], const int rank = -1)
2059: {
2060: ostringstream output;
2061: ostringstream rankStr;
2063: if (rank >= 0) {
2064: rankStr << "[" << rank << "]";
2065: }
2066: output << rankStr.str() << name << " = " << std::endl;
2067: for(int r = 0; r < rows; r++) {
2068: if (r == 0) {
2069: output << rankStr.str() << " /";
2070: } else if (r == rows-1) {
2071: output << rankStr.str() << " \\";
2072: } else {
2073: output << rankStr.str() << " |";
2074: }
2075: for(int c = 0; c < cols; c++) {
2076: output << " " << matrix[r*cols+c];
2077: }
2078: if (r == 0) {
2079: output << " \\" << std::endl;
2080: } else if (r == rows-1) {
2081: output << " /" << std::endl;
2082: } else {
2083: output << " |" << std::endl;
2084: }
2085: }
2086: return output.str();
2087: };
2088: };
2089: class MeshBuilder {
2090: typedef ALE::Mesh Mesh;
2091: public:
2094: /*
2095: Simple square boundary:
2097: 18--5-17--4--16
2098: | | |
2099: 6 10 3
2100: | | |
2101: 19-11-20--9--15
2102: | | |
2103: 7 8 2
2104: | | |
2105: 12--0-13--1--14
2106: */
2107: static Obj<ALE::Mesh> createSquareBoundary(const MPI_Comm comm, const double lower[], const double upper[], const int edges[], const int debug = 0) {
2108: Obj<Mesh> mesh = new Mesh(comm, 1, debug);
2109: int numVertices = (edges[0]+1)*(edges[1]+1);
2110: int numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
2111: double *coords = new double[numVertices*2];
2112: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2113: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2114: int order = 0;
2116: mesh->setSieve(sieve);
2117: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2118: if (mesh->commRank() == 0) {
2119: /* Create sieve and ordering */
2120: for(int v = numEdges; v < numEdges+numVertices; v++) {
2121: vertices[v-numEdges] = Mesh::point_type(v);
2122: }
2123: for(int vy = 0; vy <= edges[1]; vy++) {
2124: for(int ex = 0; ex < edges[0]; ex++) {
2125: Mesh::point_type edge(vy*edges[0] + ex);
2126: int vertex = vy*(edges[0]+1) + ex;
2128: sieve->addArrow(vertices[vertex+0], edge, order++);
2129: sieve->addArrow(vertices[vertex+1], edge, order++);
2130: if ((vy == 0) || (vy == edges[1])) {
2131: mesh->setValue(markers, edge, 1);
2132: mesh->setValue(markers, vertices[vertex], 1);
2133: if (ex == edges[0]-1) {
2134: mesh->setValue(markers, vertices[vertex+1], 1);
2135: }
2136: }
2137: }
2138: }
2139: for(int vx = 0; vx <= edges[0]; vx++) {
2140: for(int ey = 0; ey < edges[1]; ey++) {
2141: Mesh::point_type edge(vx*edges[1] + ey + edges[0]*(edges[1]+1));
2142: int vertex = ey*(edges[0]+1) + vx;
2144: sieve->addArrow(vertices[vertex], edge, order++);
2145: sieve->addArrow(vertices[vertex+edges[0]+1], edge, order++);
2146: if ((vx == 0) || (vx == edges[0])) {
2147: mesh->setValue(markers, edge, 1);
2148: mesh->setValue(markers, vertices[vertex], 1);
2149: if (ey == edges[1]-1) {
2150: mesh->setValue(markers, vertices[vertex+edges[0]+1], 1);
2151: }
2152: }
2153: }
2154: }
2155: }
2156: mesh->stratify();
2157: for(int vy = 0; vy <= edges[1]; ++vy) {
2158: for(int vx = 0; vx <= edges[0]; ++vx) {
2159: coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
2160: coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
2161: }
2162: }
2163: delete [] vertices;
2164: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2165: return mesh;
2166: };
2169: /*
2170: Simple square boundary:
2172: 18--5-17--4--16
2173: | | |
2174: 6 10 3
2175: | | |
2176: 19-11-20--9--15
2177: | | |
2178: 7 8 2
2179: | | |
2180: 12--0-13--1--14
2181: */
2182: static Obj<ALE::Mesh> createParticleInSquareBoundary(const MPI_Comm comm, const double lower[], const double upper[], const int edges[], const double radius, const int partEdges, const int debug = 0) {
2183: Obj<Mesh> mesh = new Mesh(comm, 1, debug);
2184: const int numSquareVertices = (edges[0]+1)*(edges[1]+1);
2185: const int numVertices = numSquareVertices + partEdges;
2186: const int numSquareEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1];
2187: const int numEdges = numSquareEdges + partEdges;
2188: double *coords = new double[numVertices*2];
2189: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2190: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2191: int order = 0;
2193: mesh->setSieve(sieve);
2194: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2195: if (mesh->commRank() == 0) {
2196: /* Create sieve and ordering */
2197: for(int v = numEdges; v < numEdges+numVertices; v++) {
2198: vertices[v-numEdges] = Mesh::point_type(v);
2199: }
2200: // Make square
2201: for(int vy = 0; vy <= edges[1]; vy++) {
2202: for(int ex = 0; ex < edges[0]; ex++) {
2203: Mesh::point_type edge(vy*edges[0] + ex);
2204: int vertex = vy*(edges[0]+1) + ex;
2206: sieve->addArrow(vertices[vertex+0], edge, order++);
2207: sieve->addArrow(vertices[vertex+1], edge, order++);
2208: if ((vy == 0) || (vy == edges[1])) {
2209: mesh->setValue(markers, edge, 1);
2210: mesh->setValue(markers, vertices[vertex], 1);
2211: if (ex == edges[0]-1) {
2212: mesh->setValue(markers, vertices[vertex+1], 1);
2213: }
2214: }
2215: }
2216: }
2217: for(int vx = 0; vx <= edges[0]; vx++) {
2218: for(int ey = 0; ey < edges[1]; ey++) {
2219: Mesh::point_type edge(vx*edges[1] + ey + edges[0]*(edges[1]+1));
2220: int vertex = ey*(edges[0]+1) + vx;
2222: sieve->addArrow(vertices[vertex], edge, order++);
2223: sieve->addArrow(vertices[vertex+edges[0]+1], edge, order++);
2224: if ((vx == 0) || (vx == edges[0])) {
2225: mesh->setValue(markers, edge, 1);
2226: mesh->setValue(markers, vertices[vertex], 1);
2227: if (ey == edges[1]-1) {
2228: mesh->setValue(markers, vertices[vertex+edges[0]+1], 1);
2229: }
2230: }
2231: }
2232: }
2233: // Make particle
2234: for(int ep = 0; ep < partEdges; ++ep) {
2235: Mesh::point_type edge(numSquareEdges + ep);
2236: const int vertexA = numSquareVertices + ep;
2237: const int vertexB = numSquareVertices + (ep+1)%partEdges;
2239: sieve->addArrow(vertices[vertexA], edge, order++);
2240: sieve->addArrow(vertices[vertexB], edge, order++);
2241: mesh->setValue(markers, edge, 2);
2242: mesh->setValue(markers, vertices[vertexA], 2);
2243: mesh->setValue(markers, vertices[vertexB], 2);
2244: }
2245: }
2246: mesh->stratify();
2247: for(int vy = 0; vy <= edges[1]; ++vy) {
2248: for(int vx = 0; vx <= edges[0]; ++vx) {
2249: coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx;
2250: coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy;
2251: }
2252: }
2253: const double centroidX = 0.5*(upper[0] + lower[0]);
2254: const double centroidY = 0.5*(upper[1] + lower[1]);
2255: for(int vp = 0; vp < partEdges; ++vp) {
2256: const double rad = 2.0*PETSC_PI*vp/partEdges;
2257: coords[(numSquareVertices+vp)*2+0] = centroidX + radius*cos(rad);
2258: coords[(numSquareVertices+vp)*2+1] = centroidY + radius*sin(rad);
2259: }
2260: delete [] vertices;
2261: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2262: return mesh;
2263: };
2266: /*
2267: Simple boundary with reentrant singularity:
2269: 12--5-11
2270: | |
2271: | 4
2272: | |
2273: 6 10--3--9
2274: | |
2275: | 2
2276: | |
2277: 7-----1-----8
2278: */
2279: static Obj<ALE::Mesh> createReentrantBoundary(const MPI_Comm comm, const double lower[], const double upper[], double notchpercent[], const int debug = 0) {
2280: Obj<Mesh> mesh = new Mesh(comm, 1, debug);
2281: int numVertices = 6;
2282: int numEdges = numVertices;
2283: double *coords = new double[numVertices*2];
2284: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2285: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2287: mesh->setSieve(sieve);
2288: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2289: if (mesh->commRank() == 0) {
2290: /* Create sieve and ordering */
2291: for (int b = 0; b < numVertices; b++) {
2292: sieve->addArrow(numEdges+b, b);
2293: sieve->addArrow(numEdges+b, (b+1)%numVertices);
2294: mesh->setValue(markers, b, 1);
2295: mesh->setValue(markers, b+numVertices, 1);
2296: }
2297: coords[0] = upper[0];
2298: coords[1] = lower[1];
2300: coords[2] = lower[0];
2301: coords[3] = lower[1];
2302:
2303: coords[4] = lower[0];
2304: coords[5] = notchpercent[1]*lower[1] + (1 - notchpercent[1])*upper[1];
2305:
2306: coords[6] = notchpercent[0]*upper[0] + (1 - notchpercent[0])*lower[0];
2307: coords[7] = notchpercent[1]*lower[1] + (1 - notchpercent[1])*upper[1];
2308:
2309:
2310: coords[8] = notchpercent[0]*upper[0] + (1 - notchpercent[0])*lower[0];
2311: coords[9] = upper[1];
2313: coords[10] = upper[0];
2314: coords[11] = upper[1];
2315: mesh->stratify();
2316: }
2317: delete [] vertices;
2318: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2319: return mesh;
2320: }
2324: /*
2325: Circular boundary with reentrant singularity:
2327: ---1
2328: -- |
2329: - |
2330: | |
2331: | 0-----n
2332: | |
2333: - -
2334: -- --
2335: -------
2336: */
2337: static Obj<ALE::Mesh> createCircularReentrantBoundary(const MPI_Comm comm, const int segments, const double radius, const double arc_percent, const int debug = 0) {
2338: Obj<Mesh> mesh = new Mesh(comm, 1, debug);
2339: int numVertices = segments+2;
2340: int numEdges = numVertices;
2341: double *coords = new double[numVertices*2];
2342: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2343: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2345: mesh->setSieve(sieve);
2346: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2347: if (mesh->commRank() == 0) {
2348: /* Create sieve and ordering */
2350: int startvertex = 1;
2351: if (arc_percent < 1.) {
2352: coords[0] = 0.;
2353: coords[1] = 0.;
2354: } else {
2355: numVertices = segments;
2356: numEdges = numVertices;
2357: startvertex = 0;
2358: }
2360: for (int b = 0; b < numVertices; b++) {
2361: sieve->addArrow(numEdges+b, b);
2362: sieve->addArrow(numEdges+b, (b+1)%numVertices);
2363: mesh->setValue(markers, b, 1);
2364: mesh->setValue(markers, b+numVertices, 1);
2365: }
2367: double anglestep = arc_percent*2.*3.14159265/((float)segments);
2369: for (int i = startvertex; i < numVertices; i++) {
2370: coords[2*i] = radius * sin(anglestep*(i-startvertex));
2371: coords[2*i+1] = radius*cos(anglestep*(i-startvertex));
2372: }
2373: mesh->stratify();
2374: }
2375: delete [] vertices;
2376: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2377: return mesh;
2378: };
2381: static Obj<ALE::Mesh> createAnnularBoundary(const MPI_Comm comm, const int segments, const double centers[4], const double radii[2], const int debug = 0) {
2382: Obj<Mesh> mesh = new Mesh(comm, 1, debug);
2383: int numVertices = segments*2;
2384: int numEdges = numVertices;
2385: double *coords = new double[numVertices*2];
2386: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2387: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2389: mesh->setSieve(sieve);
2390: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2391: if (mesh->commRank() == 0) {
2392: for (int e = 0; e < segments; ++e) {
2393: sieve->addArrow(numEdges+e, e);
2394: sieve->addArrow(numEdges+(e+1)%segments, e);
2395: sieve->addArrow(numEdges+segments+e, e+segments);
2396: sieve->addArrow(numEdges+segments+(e+1)%segments, e+segments);
2397: mesh->setValue(markers, e, 1);
2398: mesh->setValue(markers, e+segments, 1);
2399: mesh->setValue(markers, e+numEdges, 1);
2400: mesh->setValue(markers, e+numEdges+segments, 1);
2401: }
2402: const double anglestep = 2.0*M_PI/segments;
2404: for (int v = 0; v < segments; ++v) {
2405: coords[v*2] = centers[0] + radii[0]*cos(anglestep*v);
2406: coords[v*2+1] = centers[1] + radii[0]*sin(anglestep*v);
2407: coords[(v+segments)*2] = centers[2] + radii[1]*cos(anglestep*v);
2408: coords[(v+segments)*2+1] = centers[3] + radii[1]*sin(anglestep*v);
2409: }
2410: mesh->addHole(¢ers[2]);
2411: }
2412: mesh->stratify();
2413: delete [] vertices;
2414: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2415: return mesh;
2416: };
2419: /*
2420: Simple cubic boundary:
2422: 30----31-----32
2423: | | |
2424: | 3 | 2 |
2425: | | |
2426: 27----28-----29
2427: | | |
2428: | 0 | 1 |
2429: | | |
2430: 24----25-----26
2431: */
2432: static Obj<Mesh> createCubeBoundary(const MPI_Comm comm, const double lower[], const double upper[], const int faces[], const int debug = 0) {
2433: Obj<Mesh> mesh = new Mesh(comm, 2, debug);
2434: int numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
2435: int numFaces = 6;
2436: double *coords = new double[numVertices*3];
2437: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2438: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2439: int order = 0;
2441: mesh->setSieve(sieve);
2442: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2443: if (mesh->commRank() == 0) {
2444: /* Create sieve and ordering */
2445: for(int v = numFaces; v < numFaces+numVertices; v++) {
2446: vertices[v-numFaces] = Mesh::point_type(v);
2447: mesh->setValue(markers, vertices[v-numFaces], 1);
2448: }
2449: {
2450: // Side 0 (Front)
2451: Mesh::point_type face(0);
2452: sieve->addArrow(vertices[0], face, order++);
2453: sieve->addArrow(vertices[1], face, order++);
2454: sieve->addArrow(vertices[2], face, order++);
2455: sieve->addArrow(vertices[3], face, order++);
2456: mesh->setValue(markers, face, 1);
2457: }
2458: {
2459: // Side 1 (Back)
2460: Mesh::point_type face(1);
2461: sieve->addArrow(vertices[5], face, order++);
2462: sieve->addArrow(vertices[4], face, order++);
2463: sieve->addArrow(vertices[7], face, order++);
2464: sieve->addArrow(vertices[6], face, order++);
2465: mesh->setValue(markers, face, 1);
2466: }
2467: {
2468: // Side 2 (Bottom)
2469: Mesh::point_type face(2);
2470: sieve->addArrow(vertices[4], face, order++);
2471: sieve->addArrow(vertices[5], face, order++);
2472: sieve->addArrow(vertices[1], face, order++);
2473: sieve->addArrow(vertices[0], face, order++);
2474: mesh->setValue(markers, face, 1);
2475: }
2476: {
2477: // Side 3 (Top)
2478: Mesh::point_type face(3);
2479: sieve->addArrow(vertices[3], face, order++);
2480: sieve->addArrow(vertices[2], face, order++);
2481: sieve->addArrow(vertices[6], face, order++);
2482: sieve->addArrow(vertices[7], face, order++);
2483: mesh->setValue(markers, face, 1);
2484: }
2485: {
2486: // Side 4 (Left)
2487: Mesh::point_type face(4);
2488: sieve->addArrow(vertices[4], face, order++);
2489: sieve->addArrow(vertices[0], face, order++);
2490: sieve->addArrow(vertices[3], face, order++);
2491: sieve->addArrow(vertices[7], face, order++);
2492: mesh->setValue(markers, face, 1);
2493: }
2494: {
2495: // Side 5 (Right)
2496: Mesh::point_type face(5);
2497: sieve->addArrow(vertices[1], face, order++);
2498: sieve->addArrow(vertices[5], face, order++);
2499: sieve->addArrow(vertices[6], face, order++);
2500: sieve->addArrow(vertices[2], face, order++);
2501: mesh->setValue(markers, face, 1);
2502: }
2503: }
2504: mesh->stratify();
2505: #if 0
2506: for(int vz = 0; vz <= edges[2]; ++vz) {
2507: for(int vy = 0; vy <= edges[1]; ++vy) {
2508: for(int vx = 0; vx <= edges[0]; ++vx) {
2509: coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
2510: coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
2511: coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
2512: }
2513: }
2514: }
2515: #else
2516: coords[0*3+0] = lower[0];
2517: coords[0*3+1] = lower[1];
2518: coords[0*3+2] = upper[2];
2519: coords[1*3+0] = upper[0];
2520: coords[1*3+1] = lower[1];
2521: coords[1*3+2] = upper[2];
2522: coords[2*3+0] = upper[0];
2523: coords[2*3+1] = upper[1];
2524: coords[2*3+2] = upper[2];
2525: coords[3*3+0] = lower[0];
2526: coords[3*3+1] = upper[1];
2527: coords[3*3+2] = upper[2];
2528: coords[4*3+0] = lower[0];
2529: coords[4*3+1] = lower[1];
2530: coords[4*3+2] = lower[2];
2531: coords[5*3+0] = upper[0];
2532: coords[5*3+1] = lower[1];
2533: coords[5*3+2] = lower[2];
2534: coords[6*3+0] = upper[0];
2535: coords[6*3+1] = upper[1];
2536: coords[6*3+2] = lower[2];
2537: coords[7*3+0] = lower[0];
2538: coords[7*3+1] = upper[1];
2539: coords[7*3+2] = lower[2];
2540: #endif
2541: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2542: return mesh;
2543: };
2547: /* v0
2548: / \
2549: / \
2550: 2 / 4 \ 1
2551: / \
2552: / v12 \
2553: v6|\ /|\ /|v5
2554: | \v8 | v7/ | z
2555: | |7 |8 | | |
2556: | |v13\ | | <-v4 / \
2557: | v9/ 9 \v10| x y
2558: v1 | 5 \ / 6 |v2
2559: \ \ / /
2560: \ v11 /
2561: \ | /
2562: 3 \ | /
2563: \|/
2564: v3
2565: */
2566: static Obj<Mesh> createFicheraCornerBoundary(const MPI_Comm comm, const double lower[], const double upper[], const double offset[], const int debug = 0) {
2567: Obj<Mesh> mesh = new Mesh(comm, 2, debug);
2568: const int nVertices = 14;
2569: const int nFaces = 12;
2570: double ilower[3];
2571: ilower[0] = lower[0]*(1. - offset[0]) + upper[0]*offset[0];
2572: ilower[1] = lower[1]*(1. - offset[1]) + upper[1]*offset[1];
2573: ilower[2] = lower[2]*(1. - offset[2]) + upper[2]*offset[2];
2574: double coords[nVertices*3];
2575: //outer square-triplet
2576: coords[0*3+0] = lower[0];
2577: coords[0*3+1] = lower[1];
2578: coords[0*3+2] = upper[2];
2579: coords[1*3+0] = upper[0];
2580: coords[1*3+1] = lower[1];
2581: coords[1*3+2] = lower[2];
2582: coords[2*3+0] = lower[0];
2583: coords[2*3+1] = upper[1];
2584: coords[2*3+2] = lower[2];
2585: coords[3*3+0] = upper[0];
2586: coords[3*3+1] = upper[1];
2587: coords[3*3+2] = lower[2];
2588: coords[4*3+0] = lower[0];
2589: coords[4*3+1] = lower[1];
2590: coords[4*3+2] = lower[2];
2591: coords[5*3+0] = lower[0];
2592: coords[5*3+1] = upper[1];
2593: coords[5*3+2] = upper[2];
2594: coords[6*3+0] = upper[0];
2595: coords[6*3+1] = lower[1];
2596: coords[6*3+2] = upper[2];
2598: //inner square-triplet
2599: coords[7*3+0] = ilower[0];
2600: coords[7*3+1] = upper[1];
2601: coords[7*3+2] = upper[2];
2602: coords[8*3+0] = upper[0];
2603: coords[8*3+1] = ilower[1];
2604: coords[8*3+2] = upper[2];
2605: coords[9*3+0] = upper[0];
2606: coords[9*3+1] = ilower[1];
2607: coords[9*3+2] = ilower[2];
2608: coords[10*3+0] = ilower[0];
2609: coords[10*3+1] = upper[1];
2610: coords[10*3+2] = ilower[2];
2611: coords[11*3+0] = upper[0];
2612: coords[11*3+1] = upper[1];
2613: coords[11*3+2] = ilower[2];
2614: coords[12*3+0] = ilower[0];
2615: coords[12*3+1] = ilower[1];
2616: coords[12*3+2] = upper[2];
2617: coords[13*3+0] = ilower[0];
2618: coords[13*3+1] = ilower[1];
2619: coords[13*3+2] = ilower[2];
2621:
2622: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2623: mesh->setSieve(sieve);
2624: Mesh::point_type p[nVertices];
2625: Mesh::point_type f[nFaces];
2626: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2627: for (int i = 0; i < nVertices; i++) {
2628: p[i] = Mesh::point_type(i+nFaces);
2629: mesh->setValue(markers, p[i], 1);
2630: }
2631: for (int i = 0; i < nFaces; i++) {
2632: f[i] = Mesh::point_type(i);
2633: }
2634: int order = 0;
2635: //assemble the larger square sides
2636: sieve->addArrow(p[0], f[0], order++);
2637: sieve->addArrow(p[5], f[0], order++);
2638: sieve->addArrow(p[2], f[0], order++);
2639: sieve->addArrow(p[4], f[0], order++);
2640: mesh->setValue(markers, f[0], 1);
2642: sieve->addArrow(p[0], f[1], order++);
2643: sieve->addArrow(p[4], f[1], order++);
2644: sieve->addArrow(p[1], f[1], order++);
2645: sieve->addArrow(p[6], f[1], order++);
2646: mesh->setValue(markers, f[1], 1);
2648: sieve->addArrow(p[4], f[2], order++);
2649: sieve->addArrow(p[1], f[2], order++);
2650: sieve->addArrow(p[3], f[2], order++);
2651: sieve->addArrow(p[2], f[2], order++);
2652: mesh->setValue(markers, f[2], 1);
2653:
2654: //assemble the L-shaped sides
2656: sieve->addArrow(p[0], f[3], order++);
2657: sieve->addArrow(p[12], f[3], order++);
2658: sieve->addArrow(p[7], f[3], order++);
2659: sieve->addArrow(p[5], f[3], order++);
2660: mesh->setValue(markers, f[3], 1);
2662: sieve->addArrow(p[0], f[4], order++);
2663: sieve->addArrow(p[12],f[4], order++);
2664: sieve->addArrow(p[8], f[4], order++);
2665: sieve->addArrow(p[6], f[4], order++);
2666: mesh->setValue(markers, f[4], 1);
2668: sieve->addArrow(p[9], f[5], order++);
2669: sieve->addArrow(p[1], f[5], order++);
2670: sieve->addArrow(p[3], f[5], order++);
2671: sieve->addArrow(p[11], f[5], order++);
2672: mesh->setValue(markers, f[5], 1);
2674: sieve->addArrow(p[9], f[6], order++);
2675: sieve->addArrow(p[1], f[6], order++);
2676: sieve->addArrow(p[6], f[6], order++);
2677: sieve->addArrow(p[8], f[6], order++);
2678: mesh->setValue(markers, f[6], 1);
2680: sieve->addArrow(p[10], f[7], order++);
2681: sieve->addArrow(p[2], f[7], order++);
2682: sieve->addArrow(p[5], f[7], order++);
2683: sieve->addArrow(p[7], f[7], order++);
2684: mesh->setValue(markers, f[7], 1);
2686: sieve->addArrow(p[10], f[8], order++);
2687: sieve->addArrow(p[2], f[8], order++);
2688: sieve->addArrow(p[3], f[8], order++);
2689: sieve->addArrow(p[11], f[8], order++);
2690: mesh->setValue(markers, f[8], 1);
2692: //assemble the smaller square sides
2694: sieve->addArrow(p[13], f[9], order++);
2695: sieve->addArrow(p[10], f[9], order++);
2696: sieve->addArrow(p[11], f[9], order++);
2697: sieve->addArrow(p[9], f[9], order++);
2698: mesh->setValue(markers, f[9], 1);
2700: sieve->addArrow(p[12], f[10], order++);
2701: sieve->addArrow(p[7], f[10], order++);
2702: sieve->addArrow(p[10], f[10], order++);
2703: sieve->addArrow(p[13], f[10], order++);
2704: mesh->setValue(markers, f[10], 1);
2706: sieve->addArrow(p[8], f[11], order++);
2707: sieve->addArrow(p[12], f[11], order++);
2708: sieve->addArrow(p[13], f[11], order++);
2709: sieve->addArrow(p[9], f[11], order++);
2710: mesh->setValue(markers, f[11], 1);
2712: mesh->stratify();
2713: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
2714: Obj<Mesh::real_section_type> coordinates = mesh->getRealSection("coordinates");
2715: //coordinates->view("coordinates");
2716: return mesh;
2718: }
2722: /*
2723: //"sphere" out a cube
2725: */
2726: #if 0
2727: static Obj<Mesh> createSphereBoundary(const MPI_Comm comm, const double radius, const int refinement, const int debug = 0) {
2728: Obj<Mesh> m = new Mesh(comm, 2, debug);
2729: Obj<Mesh::sieve_type> s = new Mesh::sieve_type(comm, debug);
2730: m->setSieve(s);
2731: Mesh::point_type p = 0;
2732: int nVertices = 8+12*(refinement)+6*(refinement)*(refinement);
2733: Mesh::point_type vertices[nVertices];
2734: double coords[3*nVertices];
2735: int nCells = 6*2*(refinement+1)*(refinement+1);
2736: double delta = 2./((double)(refinement+1));
2737: Mesh::point_type cells[nCells];
2738: for (int i = 0; i < nCells; i++) {
2739: cells[i] = p;
2740: p++;
2741: }
2742: for (int i = 0; i < nVertices; i++) {
2743: vertices[i] = p;
2744: p++;
2745: }
2746: //set up the corners;
2747: //lll
2748: coords[0*3+0] = -1.;
2749: coords[0*3+1] = -1.;
2750: coords[0*3+2] = -1.;
2751: //llh
2752: coords[1*3+0] = -1.;
2753: coords[1*3+1] = -1.;
2754: coords[1*3+2] = 1.;
2755: //lhh
2756: coords[2*3+0] = -1.;
2757: coords[2*3+1] = 1.;
2758: coords[2*3+2] = 1.;
2759: //lhl
2760: coords[3*3+0] = -1.;
2761: coords[3*3+1] = 1.;
2762: coords[3*3+2] = -1.;
2763: //hhl
2764: coords[4*3+0] = 1.;
2765: coords[4*3+1] = 1.;
2766: coords[4*3+2] = -1.;
2767: //hhh
2768: coords[5*3+0] = 1.;
2769: coords[5*3+1] = 1.;
2770: coords[5*3+2] = 1.;
2771: //hlh
2772: coords[6*3+0] = 1.;
2773: coords[6*3+1] = -1.;
2774: coords[6*3+2] = 1.;
2775: //hll
2776: coords[7*3+0] = 1.;
2777: coords[7*3+1] = -1.;
2778: coords[7*3+2] = -1.;
2779: //set up the edges (always go low to high)
2780: //xll
2781: for (int i = 0; i < refinement; i++) {
2782: coords[3*(8+0*refinement+i)+0] = -1. + delta*i;
2783: coords[3*(8+0*refinement+i)+1] = -1.;
2784: coords[3*(8+0*refinement+i)+2] = -1.;
2785: }
2786: //xlh
2787: for (int i = 0; i < refinement; i++) {
2788: coords[3*(8+1*refinement+i)+0] = -1. + delta*i;
2789: coords[3*(8+1*refinement+i)+1] = -1.;
2790: coords[3*(8+1*refinement+i)+2] = 1.;
2791: }
2792: //xhh
2793: for (int i = 0; i < refinement; i++) {
2794: coords[3*(8+2*refinement+i)+0] = -1. + delta*i;
2795: coords[3*(8+2*refinement+i)+1] = 1.;
2796: coords[3*(8+2*refinement+i)+2] = 1.;
2797: }
2798: //xhl
2799: for (int i = 0; i < refinement; i++) {
2800: coords[3*(8+3*refinement+i)+0] = -1. + delta*i;
2801: coords[3*(8+3*refinement+i)+1] = 1.;
2802: coords[3*(8+3*refinement+i)+2] = -1.;
2803: }
2804: //lxl
2805: for (int i = 0; i < refinement; i++) {
2806: coords[3*(8+4*refinement+i)+0] = -1.;
2807: coords[3*(8+4*refinement+i)+1] = -1. + delta*i;
2808: coords[3*(8+4*refinement+i)+2] = -1.;
2809: }
2810: //lxh
2811: for (int i = 0; i < refinement; i++) {
2812: coords[3*(8+5*refinement+i)+0] = -1.;
2813: coords[3*(8+5*refinement+i)+1] = -1. + delta*i;
2814: coords[3*(8+5*refinement+i)+2] = 1.;
2815: }
2816: //hxh
2817: for (int i = 0; i < refinement; i++) {
2818: coords[3*(8+6*refinement+i)+0] = 1.;
2819: coords[3*(8+6*refinement+i)+1] = -1. + delta*i;
2820: coords[3*(8+6*refinement+i)+2] = 1.;
2821: }
2822: //hxl
2823: for (int i = 0; i < refinement; i++) {
2824: coords[3*(8+7*refinement+i)+0] = 1.;
2825: coords[3*(8+7*refinement+i)+1] = -1. + delta*i;
2826: coords[3*(8+7*refinement+i)+2] = -1.;
2827: }
2828: //llx
2829: for (int i = 0; i < refinement; i++) {
2830: coords[3*(8+8*refinement+i)+0] = -1.;
2831: coords[3*(8+8*refinement+i)+1] = -1.;
2832: coords[3*(8+8*refinement+i)+2] = -1. + delta*i;
2833: }
2834: //lhx
2835: for (int i = 0; i < refinement; i++) {
2836: coords[3*(8+9*refinement+i)+0] = -1.;
2837: coords[3*(8+9*refinement+i)+1] = 1.;
2838: coords[3*(8+9*refinement+i)+2] = -1. + delta*i;
2839: }
2840: //hhx
2841: for (int i = 0; i < refinement; i++) {
2842: coords[3*(8+10*refinement+i)+0] = 1.;
2843: coords[3*(8+10*refinement+i)+1] = 1.;
2844: coords[3*(8+10*refinement+i)+2] = -1. + delta*i;
2845: }
2846: //hlx
2847: for (int i = 0; i < refinement; i++) {
2848: coords[3*(8+11*refinement+i)+0] = 1.;
2849: coords[3*(8+11*refinement+i)+1] = -1.;
2850: coords[3*(8+11*refinement+i)+2] = -1. + delta*i;
2851: }
2852: //set up the faces
2853: //lxx
2854: for (int i = 0; i < refinement; i++) for (int j = 0; j < refinement; j++) {
2855: coords[3*(8+12*refinement+0*refinement*refinement+i*refinement+j)+0] = -1.;
2856: coords[3*(8+12*refinement+0*refinement*refinement+i*refinement+j)+1] = -1. + delta*j;
2857: coords[3*(8+12*refinement+0*refinement*refinement+i*refinement+j)+2] = -1. + delta*i;
2858: }
2859: //hxx
2860: for (int i = 0; i < refinement; i++) for (int j = 0; j < refinement; j++) {
2861: coords[3*(8+12*refinement+1*refinement*refinement+i*refinement+j)+0] = 1.;
2862: coords[3*(8+12*refinement+1*refinement*refinement+i*refinement+j)+1] = -1. + delta*j;
2863: coords[3*(8+12*refinement+1*refinement*refinement+i*refinement+j)+2] = -1. + delta*i;
2864: }
2865: //xlx
2866: for (int i = 0; i < refinement; i++) for (int j = 0; j < refinement; j++) {
2867: coords[3*(8+12*refinement+2*refinement*refinement+i*refinement+j)+0] = -1. + delta*j;
2868: coords[3*(8+12*refinement+2*refinement*refinement+i*refinement+j)+1] = -1.;
2869: coords[3*(8+12*refinement+2*refinement*refinement+i*refinement+j)+2] = -1. + delta*i;
2870: }
2871: //xhx
2872: for (int i = 0; i < refinement; i++) for (int j = 0; j < refinement; j++) {
2873: coords[3*(8+12*refinement+3*refinement*refinement+i*refinement+j)+0] = -1. + delta*j;
2874: coords[3*(8+12*refinement+3*refinement*refinement+i*refinement+j)+1] = 1.;
2875: coords[3*(8+12*refinement+3*refinement*refinement+i*refinement+j)+2] = -1. + delta*i;
2876: }
2877: //xxl
2878: for (int i = 0; i < refinement; i++) for (int j = 0; j < refinement; j++) {
2879: coords[3*(8+12*refinement+4*refinement*refinement+i*refinement+j)+0] = -1.;
2880: coords[3*(8+12*refinement+4*refinement*refinement+i*refinement+j)+1] = -1. + delta*j;
2881: coords[3*(8+12*refinement+4*refinement*refinement+i*refinement+j)+2] = -1. + delta*i;
2882: }
2883: //xxh
2884: for (int i = 0; i < refinement; i++) for (int j = 0; j < refinement; j++) {
2885: coords[3*(8+12*refinement+5*refinement*refinement+i*refinement+j)+0] = 1.;
2886: coords[3*(8+12*refinement+5*refinement*refinement+i*refinement+j)+1] = -1. + delta*j;
2887: coords[3*(8+12*refinement+5*refinement*refinement+i*refinement+j)+2] = -1. + delta*i;
2888: }
2889: //stitch the corners up with the edges and the faces
2890:
2891: //stitch the edges to the faces
2892: //fill in the faces
2893: int face_offset = 8 + 12*refinement;
2894: for (int i = 0; i < 6; i++) for (int j = 0; j < refinement; j++) for (int k = 0; k < refinement; k++) {
2895: //build each square doublet
2896: }
2897: }
2899: #endif
2903: /*
2904: Simple cubic boundary:
2906: 30----31-----32
2907: | | |
2908: | 3 | 2 |
2909: | | |
2910: 27----28-----29
2911: | | |
2912: | 0 | 1 |
2913: | | |
2914: 24----25-----26
2915: */
2916: static Obj<Mesh> createParticleInCubeBoundary(const MPI_Comm comm, const double lower[], const double upper[], const int faces[], const double radius, const int thetaEdges, const int phiSlices, const int debug = 0) {
2917: Obj<Mesh> mesh = new Mesh(comm, 2, debug);
2918: const int numCubeVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1);
2919: const int numPartVertices = (thetaEdges - 1)*phiSlices + 2;
2920: const int numVertices = numCubeVertices + numPartVertices;
2921: const int numCubeFaces = 6;
2922: const int numFaces = numCubeFaces + thetaEdges*phiSlices;
2923: double *coords = new double[numVertices*3];
2924: const Obj<Mesh::sieve_type> sieve = new Mesh::sieve_type(mesh->comm(), mesh->debug());
2925: Mesh::point_type *vertices = new Mesh::point_type[numVertices];
2926: int order = 0;
2928: mesh->setSieve(sieve);
2929: const Obj<Mesh::label_type>& markers = mesh->createLabel("marker");
2930: if (mesh->commRank() == 0) {
2931: // Make cube
2932: for(int v = numFaces; v < numFaces+numVertices; v++) {
2933: vertices[v-numFaces] = Mesh::point_type(v);
2934: mesh->setValue(markers, vertices[v-numFaces], 1);
2935: }
2936: {
2937: // Side 0 (Front)
2938: Mesh::point_type face(0);
2939: sieve->addArrow(vertices[0], face, order++);
2940: sieve->addArrow(vertices[1], face, order++);
2941: sieve->addArrow(vertices[2], face, order++);
2942: sieve->addArrow(vertices[3], face, order++);
2943: mesh->setValue(markers, face, 1);
2944: }
2945: {
2946: // Side 1 (Back)
2947: Mesh::point_type face(1);
2948: sieve->addArrow(vertices[5], face, order++);
2949: sieve->addArrow(vertices[4], face, order++);
2950: sieve->addArrow(vertices[7], face, order++);
2951: sieve->addArrow(vertices[6], face, order++);
2952: mesh->setValue(markers, face, 1);
2953: }
2954: {
2955: // Side 2 (Bottom)
2956: Mesh::point_type face(2);
2957: sieve->addArrow(vertices[4], face, order++);
2958: sieve->addArrow(vertices[5], face, order++);
2959: sieve->addArrow(vertices[1], face, order++);
2960: sieve->addArrow(vertices[0], face, order++);
2961: mesh->setValue(markers, face, 1);
2962: }
2963: {
2964: // Side 3 (Top)
2965: Mesh::point_type face(3);
2966: sieve->addArrow(vertices[3], face, order++);
2967: sieve->addArrow(vertices[2], face, order++);
2968: sieve->addArrow(vertices[6], face, order++);
2969: sieve->addArrow(vertices[7], face, order++);
2970: mesh->setValue(markers, face, 1);
2971: }
2972: {
2973: // Side 4 (Left)
2974: Mesh::point_type face(4);
2975: sieve->addArrow(vertices[4], face, order++);
2976: sieve->addArrow(vertices[0], face, order++);
2977: sieve->addArrow(vertices[3], face, order++);
2978: sieve->addArrow(vertices[7], face, order++);
2979: mesh->setValue(markers, face, 1);
2980: }
2981: {
2982: // Side 5 (Right)
2983: Mesh::point_type face(5);
2984: sieve->addArrow(vertices[1], face, order++);
2985: sieve->addArrow(vertices[5], face, order++);
2986: sieve->addArrow(vertices[6], face, order++);
2987: sieve->addArrow(vertices[2], face, order++);
2988: mesh->setValue(markers, face, 1);
2989: }
2990: // Make particle
2991: for(int s = 0; s < phiSlices; ++s) {
2992: for(int ep = 0; ep < thetaEdges; ++ep) {
2993: // Vertices on each slice are 0..thetaEdges
2994: Mesh::point_type face(numCubeFaces + s*thetaEdges + ep);
2995: int vertexA = numCubeVertices + ep + 0 + s*(thetaEdges+1);
2996: int vertexB = numCubeVertices + ep + 1 + s*(thetaEdges+1);
2997: int vertexC = numCubeVertices + (ep + 1 + (s+1)*(thetaEdges+1))%((thetaEdges+1)*phiSlices);
2998: int vertexD = numCubeVertices + (ep + 0 + (s+1)*(thetaEdges+1))%((thetaEdges+1)*phiSlices);
2999: const int correction1 = (s > 0)*((s-1)*2 + 1);
3000: const int correction2 = (s < phiSlices-1)*(s*2 + 1);
3002: if ((vertexA - numCubeVertices)%(thetaEdges+1) == 0) {
3003: vertexA = vertexD = numCubeVertices;
3004: vertexB -= correction1;
3005: vertexC -= correction2;
3006: } else if ((vertexB - numCubeVertices)%(thetaEdges+1) == thetaEdges) {
3007: vertexA -= correction1;
3008: vertexD -= correction2;
3009: vertexB = vertexC = numCubeVertices + thetaEdges;
3010: } else {
3011: vertexA -= correction1;
3012: vertexB -= correction1;
3013: vertexC -= correction2;
3014: vertexD -= correction2;
3015: }
3016: if ((vertexA >= numVertices) || (vertexB >= numVertices) || (vertexC >= numVertices) || (vertexD >= numVertices)) {
3017: throw ALE::Exception("Bad vertex");
3018: }
3019: sieve->addArrow(vertices[vertexA], face, order++);
3020: sieve->addArrow(vertices[vertexB], face, order++);
3021: if (vertexB != vertexC) sieve->addArrow(vertices[vertexC], face, order++);
3022: if (vertexA != vertexD) sieve->addArrow(vertices[vertexD], face, order++);
3023: mesh->setValue(markers, face, 2);
3024: mesh->setValue(markers, vertices[vertexA], 2);
3025: mesh->setValue(markers, vertices[vertexB], 2);
3026: if (vertexB != vertexC) mesh->setValue(markers, vertices[vertexC], 2);
3027: if (vertexA != vertexD) mesh->setValue(markers, vertices[vertexD], 2);
3028: }
3029: }
3030: }
3031: mesh->stratify();
3032: #if 0
3033: for(int vz = 0; vz <= edges[2]; ++vz) {
3034: for(int vy = 0; vy <= edges[1]; ++vy) {
3035: for(int vx = 0; vx <= edges[0]; ++vx) {
3036: coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx;
3037: coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy;
3038: coords[((vz*(edges[1]+1)+vy)*(edges[0]+1)+vx)*2+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz;
3039: }
3040: }
3041: }
3042: #else
3043: coords[0*3+0] = lower[0];
3044: coords[0*3+1] = lower[1];
3045: coords[0*3+2] = upper[2];
3046: coords[1*3+0] = upper[0];
3047: coords[1*3+1] = lower[1];
3048: coords[1*3+2] = upper[2];
3049: coords[2*3+0] = upper[0];
3050: coords[2*3+1] = upper[1];
3051: coords[2*3+2] = upper[2];
3052: coords[3*3+0] = lower[0];
3053: coords[3*3+1] = upper[1];
3054: coords[3*3+2] = upper[2];
3055: coords[4*3+0] = lower[0];
3056: coords[4*3+1] = lower[1];
3057: coords[4*3+2] = lower[2];
3058: coords[5*3+0] = upper[0];
3059: coords[5*3+1] = lower[1];
3060: coords[5*3+2] = lower[2];
3061: coords[6*3+0] = upper[0];
3062: coords[6*3+1] = upper[1];
3063: coords[6*3+2] = lower[2];
3064: coords[7*3+0] = lower[0];
3065: coords[7*3+1] = upper[1];
3066: coords[7*3+2] = lower[2];
3067: #endif
3068: const double centroidX = 0.5*(upper[0] + lower[0]);
3069: const double centroidY = 0.5*(upper[1] + lower[1]);
3070: const double centroidZ = 0.5*(upper[2] + lower[2]);
3071: for(int s = 0; s < phiSlices; ++s) {
3072: for(int v = 0; v <= thetaEdges; ++v) {
3073: int vertex = numCubeVertices + v + s*(thetaEdges+1);
3074: const double theta = v*(PETSC_PI/thetaEdges);
3075: const double phi = s*(2.0*PETSC_PI/phiSlices);
3076: const int correction = (s > 0)*((s-1)*2 + 1);
3078: if ((vertex- numCubeVertices)%(thetaEdges+1) == 0) {
3079: vertex = numCubeVertices;
3080: } else if ((vertex - numCubeVertices)%(thetaEdges+1) == thetaEdges) {
3081: vertex = numCubeVertices + thetaEdges;
3082: } else {
3083: vertex -= correction;
3084: }
3085: coords[vertex*3+0] = centroidX + radius*sin(theta)*cos(phi);
3086: coords[vertex*3+1] = centroidY + radius*sin(theta)*sin(phi);
3087: coords[vertex*3+2] = centroidZ + radius*cos(theta);
3088: }
3089: }
3090: delete [] vertices;
3091: ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, mesh->getDimension()+1, coords);
3092: return mesh;
3093: };
3094: };
3095: } // namespace ALE
3097: namespace ALECompat {
3098: using ALE::Obj;
3099: template<typename Topology_>
3100: class Bundle : public ALE::ParallelObject {
3101: public:
3102: typedef Topology_ topology_type;
3103: typedef typename topology_type::patch_type patch_type;
3104: typedef typename topology_type::point_type point_type;
3105: typedef typename topology_type::sieve_type sieve_type;
3106: typedef typename sieve_type::coneArray coneArray;
3107: typedef ALECompat::New::Section<topology_type, double> real_section_type;
3108: typedef ALECompat::New::Section<topology_type, int> int_section_type;
3109: typedef ALE::MinimalArrow<point_type, point_type> arrow_type;
3110: typedef ALE::UniformSection<arrow_type, int> arrow_section_type;
3111: typedef struct {double x, y, z;} split_value;
3112: typedef ALE::pair<point_type, split_value> pair_type;
3113: typedef ALECompat::New::Section<topology_type, pair_type> pair_section_type;
3114: typedef std::map<std::string, Obj<real_section_type> > real_sections_type;
3115: typedef std::map<std::string, Obj<int_section_type> > int_sections_type;
3116: typedef std::map<std::string, Obj<pair_section_type> > pair_sections_type;
3117: typedef std::map<std::string, Obj<arrow_section_type> > arrow_sections_type;
3118: typedef typename topology_type::send_overlap_type send_overlap_type;
3119: typedef typename topology_type::recv_overlap_type recv_overlap_type;
3120: typedef typename ALECompat::New::SectionCompletion<topology_type, point_type>::topology_type comp_topology_type;
3121: typedef typename ALECompat::New::OverlapValues<send_overlap_type, comp_topology_type, point_type> send_section_type;
3122: typedef typename ALECompat::New::OverlapValues<recv_overlap_type, comp_topology_type, point_type> recv_section_type;
3123: protected:
3124: Obj<topology_type> _topology;
3125: bool _distributed;
3126: real_sections_type _realSections;
3127: int_sections_type _intSections;
3128: pair_sections_type _pairSections;
3129: arrow_sections_type _arrowSections;
3130: public:
3131: Bundle(MPI_Comm comm, int debug = 0) : ALE::ParallelObject(comm, debug), _distributed(false) {
3132: this->_topology = new topology_type(comm, debug);
3133: };
3134: Bundle(const Obj<topology_type>& topology) : ALE::ParallelObject(topology->comm(), topology->debug()), _topology(topology), _distributed(false) {};
3135: public: // Accessors
3136: bool getDistributed() const {return this->_distributed;};
3137: void setDistributed(const bool distributed) {this->_distributed = distributed;};
3138: const Obj<topology_type>& getTopology() const {return this->_topology;};
3139: void setTopology(const Obj<topology_type>& topology) {this->_topology = topology;};
3140: public:
3141: bool hasRealSection(const std::string& name) {
3142: return this->_realSections.find(name) != this->_realSections.end();
3143: };
3144: const Obj<real_section_type>& getRealSection(const std::string& name) {
3145: if (this->_realSections.find(name) == this->_realSections.end()) {
3146: Obj<real_section_type> section = new real_section_type(this->_topology);
3148: if (this->_debug) {std::cout << "Creating new real section: " << name << std::endl;}
3149: this->_realSections[name] = section;
3150: }
3151: return this->_realSections[name];
3152: };
3153: void setRealSection(const std::string& name, const Obj<real_section_type>& section) {
3154: this->_realSections[name] = section;
3155: };
3156: Obj<std::set<std::string> > getRealSections() const {
3157: Obj<std::set<std::string> > names = std::set<std::string>();
3159: for(typename real_sections_type::const_iterator s_iter = this->_realSections.begin(); s_iter != this->_realSections.end(); ++s_iter) {
3160: names->insert(s_iter->first);
3161: }
3162: return names;
3163: };
3164: bool hasIntSection(const std::string& name) {
3165: return this->_intSections.find(name) != this->_intSections.end();
3166: };
3167: const Obj<int_section_type>& getIntSection(const std::string& name) {
3168: if (this->_intSections.find(name) == this->_intSections.end()) {
3169: Obj<int_section_type> section = new int_section_type(this->_topology);
3171: if (this->_debug) {std::cout << "Creating new int section: " << name << std::endl;}
3172: this->_intSections[name] = section;
3173: }
3174: return this->_intSections[name];
3175: };
3176: void setIntSection(const std::string& name, const Obj<int_section_type>& section) {
3177: this->_intSections[name] = section;
3178: };
3179: Obj<std::set<std::string> > getIntSections() const {
3180: Obj<std::set<std::string> > names = std::set<std::string>();
3182: for(typename int_sections_type::const_iterator s_iter = this->_intSections.begin(); s_iter != this->_intSections.end(); ++s_iter) {
3183: names->insert(s_iter->first);
3184: }
3185: return names;
3186: };
3187: bool hasPairSection(const std::string& name) {
3188: return this->_pairSections.find(name) != this->_pairSections.end();
3189: };
3190: const Obj<pair_section_type>& getPairSection(const std::string& name) {
3191: if (this->_pairSections.find(name) == this->_pairSections.end()) {
3192: Obj<pair_section_type> section = new pair_section_type(this->_topology);
3194: if (this->_debug) {std::cout << "Creating new pair section: " << name << std::endl;}
3195: this->_pairSections[name] = section;
3196: }
3197: return this->_pairSections[name];
3198: };
3199: void setPairSection(const std::string& name, const Obj<pair_section_type>& section) {
3200: this->_pairSections[name] = section;
3201: };
3202: Obj<std::set<std::string> > getPairSections() const {
3203: Obj<std::set<std::string> > names = std::set<std::string>();
3205: for(typename pair_sections_type::const_iterator s_iter = this->_pairSections.begin(); s_iter != this->_pairSections.end(); ++s_iter) {
3206: names->insert(s_iter->first);
3207: }
3208: return names;
3209: };
3210: bool hasArrowSection(const std::string& name) {
3211: return this->_arrowSections.find(name) != this->_arrowSections.end();
3212: };
3213: const Obj<arrow_section_type>& getArrowSection(const std::string& name) {
3214: if (this->_arrowSections.find(name) == this->_arrowSections.end()) {
3215: Obj<arrow_section_type> section = new arrow_section_type(this->comm(), this->debug());
3217: if (this->_debug) {std::cout << "Creating new arrow section: " << name << std::endl;}
3218: this->_arrowSections[name] = section;
3219: }
3220: return this->_arrowSections[name];
3221: };
3222: void setArrowSection(const std::string& name, const Obj<arrow_section_type>& section) {
3223: this->_arrowSections[name] = section;
3224: };
3225: Obj<std::set<std::string> > getArrowSections() const {
3226: Obj<std::set<std::string> > names = std::set<std::string>();
3228: for(typename arrow_sections_type::const_iterator s_iter = this->_arrowSections.begin(); s_iter != this->_arrowSections.end(); ++s_iter) {
3229: names->insert(s_iter->first);
3230: }
3231: return names;
3232: };
3233: public: // Printing
3234: friend std::ostream& operator<<(std::ostream& os, const split_value& s) {
3235: os << "(" << s.x << ", "<< s.y << ", "<< s.z << ")";
3236: return os;
3237: };
3238: public: // Adapter
3239: const Obj<sieve_type>& getSieve() {return this->_topology->getPatch(0);};
3240: int height() {return 2;};
3241: int depth() {return 2;};
3242: };
3244: class Discretization : public ALE::ParallelObject {
3245: protected:
3246: std::map<int,int> _dim2dof;
3247: std::map<int,int> _dim2class;
3248: public:
3249: Discretization(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
3250: ~Discretization() {};
3251: public:
3252: const double *getQuadraturePoints() {return NULL;};
3253: const double *getQuadratureWeights() {return NULL;};
3254: const double *getBasis() {return NULL;};
3255: const double *getBasisDerivatives() {return NULL;};
3256: int getNumDof(const int dim) {return this->_dim2dof[dim];};
3257: void setNumDof(const int dim, const int numDof) {this->_dim2dof[dim] = numDof;};
3258: int getDofClass(const int dim) {return this->_dim2class[dim];};
3259: void setDofClass(const int dim, const int dofClass) {this->_dim2class[dim] = dofClass;};
3260: };
3262: class BoundaryCondition : public ALE::ParallelObject {
3263: protected:
3264: std::string _labelName;
3265: double (*_func)(const double []);
3266: public:
3267: BoundaryCondition(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
3268: ~BoundaryCondition() {};
3269: public:
3270: const std::string& getLabelName() {return this->_labelName;};
3271: void setLabelName(const std::string& name) {this->_labelName = name;};
3272: void setFunction(double (*func)(const double [])) {this->_func = func;};
3273: public:
3274: double evaluate(const double coords[]) {return this->_func(coords);};
3275: };
3277: class Mesh : public Bundle<ALE::Topology<int, ALE::Sieve<int,int,int> > > {
3278: public:
3279: typedef int point_type;
3280: typedef ALE::Sieve<point_type,int,int> sieve_type;
3281: typedef ALE::Topology<int, sieve_type> topology_type;
3282: typedef topology_type::patch_type patch_type;
3283: typedef Bundle<topology_type> base_type;
3284: typedef ALECompat::New::NumberingFactory<topology_type> NumberingFactory;
3285: typedef NumberingFactory::numbering_type numbering_type;
3286: typedef NumberingFactory::order_type order_type;
3287: typedef base_type::send_overlap_type send_overlap_type;
3288: typedef base_type::recv_overlap_type recv_overlap_type;
3289: typedef base_type::send_section_type send_section_type;
3290: typedef base_type::recv_section_type recv_section_type;
3291: typedef base_type::real_sections_type real_sections_type;
3292: // PCICE BC
3293: typedef struct {double rho,u,v,p;} bc_value_type;
3294: typedef std::map<int, bc_value_type> bc_values_type;
3295: // PyLith BC
3296: typedef ALECompat::New::Section<topology_type, ALE::pair<int,double> > foliated_section_type;
3297: protected:
3298: int _dim;
3299: Obj<NumberingFactory> _factory;
3300: // PCICE BC
3301: bc_values_type _bcValues;
3302: // PyLith BC
3303: Obj<foliated_section_type> _boundaries;
3304: // Discretization
3305: Obj<Discretization> _discretization;
3306: Obj<BoundaryCondition> _boundaryCondition;
3307: public:
3308: Mesh(MPI_Comm comm, int dim, int debug = 0) : Bundle<ALE::Topology<int, ALE::Sieve<int,int,int> > >(comm, debug), _dim(dim) {
3309: this->_factory = NumberingFactory::singleton(debug);
3310: this->_boundaries = NULL;
3311: this->_discretization = new Discretization(comm, debug);
3312: this->_boundaryCondition = new BoundaryCondition(comm, debug);
3313: };
3314: Mesh(const Obj<topology_type>& topology, int dim) : Bundle<ALE::Topology<int, ALE::Sieve<int,int,int> > >(topology), _dim(dim) {
3315: this->_factory = NumberingFactory::singleton(topology->debug());
3316: this->_boundaries = NULL;
3317: };
3318: public: // Accessors
3319: int getDimension() const {return this->_dim;};
3320: void setDimension(const int dim) {this->_dim = dim;};
3321: const Obj<NumberingFactory>& getFactory() {return this->_factory;};
3322: const Obj<Discretization>& getDiscretization() {return this->_discretization;};
3323: void setDiscretization(const Obj<Discretization>& discretization) {this->_discretization = discretization;};
3324: const Obj<BoundaryCondition>& getBoundaryCondition() {return this->_boundaryCondition;};
3325: void setBoundaryCondition(const Obj<BoundaryCondition>& boundaryCondition) {this->_boundaryCondition = boundaryCondition;};
3326: public: // Mesh geometry
3327: #if 0
3328: void computeTriangleGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
3329: const double *coords = this->restrict(coordinates, e);
3330: const int dim = 2;
3331: double invDet;
3333: if (v0) {
3334: for(int d = 0; d < dim; d++) {
3335: v0[d] = coords[d];
3336: }
3337: }
3338: if (J) {
3339: for(int d = 0; d < dim; d++) {
3340: for(int f = 0; f < dim; f++) {
3341: J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
3342: }
3343: }
3344: detJ = J[0]*J[3] - J[1]*J[2];
3345: }
3346: if (invJ) {
3347: invDet = 1.0/detJ;
3348: invJ[0] = invDet*J[3];
3349: invJ[1] = -invDet*J[1];
3350: invJ[2] = -invDet*J[2];
3351: invJ[3] = invDet*J[0];
3352: }
3353: };
3354: static void computeTetrahedronGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
3355: const patch_type patch = 0;
3356: const double *coords = coordinates->restrict(patch, e);
3357: const int dim = 3;
3358: double invDet;
3360: if (v0) {
3361: for(int d = 0; d < dim; d++) {
3362: v0[d] = coords[d];
3363: }
3364: }
3365: if (J) {
3366: for(int d = 0; d < dim; d++) {
3367: for(int f = 0; f < dim; f++) {
3368: J[d*dim+f] = 0.5*(coords[(f+1)*dim+d] - coords[0*dim+d]);
3369: }
3370: }
3371: detJ = J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
3372: J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
3373: J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
3374: }
3375: if (invJ) {
3376: invDet = 1.0/detJ;
3377: invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
3378: invJ[0*3+1] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
3379: invJ[0*3+2] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
3380: invJ[1*3+0] = invDet*(J[0*3+1]*J[2*3+2] - J[0*3+2]*J[2*3+1]);
3381: invJ[1*3+1] = invDet*(J[0*3+2]*J[2*3+0] - J[0*3+0]*J[2*3+2]);
3382: invJ[1*3+2] = invDet*(J[0*3+0]*J[2*3+1] - J[0*3+1]*J[2*3+0]);
3383: invJ[2*3+0] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
3384: invJ[2*3+1] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
3385: invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
3386: }
3387: };
3388: void computeElementGeometry(const Obj<real_section_type>& coordinates, const point_type& e, double v0[], double J[], double invJ[], double& detJ) {
3389: if (this->_dim == 2) {
3390: computeTriangleGeometry(coordinates, e, v0, J, invJ, detJ);
3391: } else if (this->_dim == 3) {
3392: computeTetrahedronGeometry(coordinates, e, v0, J, invJ, detJ);
3393: } else {
3394: throw ALE::Exception("Unsupport dimension for element geometry computation");
3395: }
3396: };
3397: double getMaxVolume() {
3398: const topology_type::sheaf_type& patches = this->getTopology()->getPatches();
3399: double v0[3], J[9], invJ[9], detJ, maxVolume = 0.0;
3401: for(topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
3402: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
3403: const Obj<topology_type::label_sequence>& cells = this->getTopology()->heightStratum(p_iter->first, 0);
3405: for(topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
3406: this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
3407: maxVolume = std::max(maxVolume, detJ);
3408: }
3409: }
3410: return maxVolume;
3411: };
3412: // Find the cell in which this point lies (stupid algorithm)
3413: point_type locatePoint_2D(const patch_type& patch, const real_section_type::value_type point[]) {
3414: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
3415: const Obj<topology_type::label_sequence>& cells = this->getTopology()->heightStratum(patch, 0);
3416: const int embedDim = 2;
3417: double v0[2], J[4], invJ[4], detJ;
3419: for(topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
3420: this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
3421: double xi = invJ[0*embedDim+0]*(point[0] - v0[0]) + invJ[0*embedDim+1]*(point[1] - v0[1]);
3422: double eta = invJ[1*embedDim+0]*(point[0] - v0[0]) + invJ[1*embedDim+1]*(point[1] - v0[1]);
3424: if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 1.0)) {
3425: return *c_iter;
3426: }
3427: }
3428: throw ALE::Exception("Could not locate point");
3429: };
3430: // Assume a simplex and 3D
3431: point_type locatePoint_3D(const patch_type& patch, const real_section_type::value_type point[]) {
3432: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
3433: const Obj<topology_type::label_sequence>& cells = this->getTopology()->heightStratum(patch, 0);
3434: const int embedDim = 3;
3435: double v0[3], J[9], invJ[9], detJ;
3437: for(topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
3438: this->computeElementGeometry(coordinates, *c_iter, v0, J, invJ, detJ);
3439: 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]);
3440: 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]);
3441: 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]);
3443: if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 1.0)) {
3444: return *c_iter;
3445: }
3446: }
3447: throw ALE::Exception("Could not locate point");
3448: };
3449: point_type locatePoint(const patch_type& patch, const real_section_type::value_type point[]) {
3450: if (this->_dim == 2) {
3451: return locatePoint_2D(patch, point);
3452: } else if (this->_dim == 3) {
3453: return locatePoint_3D(patch, point);
3454: } else {
3455: throw ALE::Exception("No point location for mesh dimension");
3456: }
3457: };
3458: #endif
3459: // Only works in 2D
3460: int orientation(const patch_type& patch, const point_type& cell) {
3461: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
3462: const Obj<topology_type>& topology = this->getTopology();
3463: const Obj<sieve_type>& sieve = topology->getPatch(patch);
3464: const Obj<sieve_type::coneArray>& cone = sieve->nCone(cell, topology->depth());
3465: sieve_type::coneArray::iterator cBegin = cone->begin();
3466: real_section_type::value_type root[2];
3467: real_section_type::value_type vA[2];
3468: real_section_type::value_type vB[2];
3470: const real_section_type::value_type *coords = coordinates->restrictPoint(patch, *cBegin);
3471: root[0] = coords[0];
3472: root[1] = coords[1];
3473: ++cBegin;
3474: coords = coordinates->restrictPoint(patch, *cBegin);
3475: vA[0] = coords[0] - root[0];
3476: vA[1] = coords[1] - root[1];
3477: ++cBegin;
3478: coords = coordinates->restrictPoint(patch, *cBegin);
3479: vB[0] = coords[0] - root[0];
3480: vB[1] = coords[1] - root[1];
3481: double det = vA[0]*vB[1] - vA[1]*vB[0];
3482: if (det > 0.0) return 1;
3483: if (det < 0.0) return -1;
3484: return 0;
3485: };
3486: public: // BC values for PCICE
3487: const bc_value_type& getBCValue(const int bcFunc) {
3488: return this->_bcValues[bcFunc];
3489: };
3490: void setBCValue(const int bcFunc, const bc_value_type& value) {
3491: this->_bcValues[bcFunc] = value;
3492: };
3493: bc_values_type& getBCValues() {
3494: return this->_bcValues;
3495: };
3496: void distributeBCValues() {
3497: int size = this->_bcValues.size();
3499: MPI_Bcast(&size, 1, MPI_INT, 0, this->comm());
3500: if (this->commRank()) {
3501: for(int bc = 0; bc < size; ++bc) {
3502: int funcNum;
3503: bc_value_type funcVal;
3505: MPI_Bcast((void *) &funcNum, 1, MPI_INT, 0, this->comm());
3506: MPI_Bcast((void *) &funcVal, 4, MPI_DOUBLE, 0, this->comm());
3507: this->_bcValues[funcNum] = funcVal;
3508: }
3509: } else {
3510: for(bc_values_type::iterator bc_iter = this->_bcValues.begin(); bc_iter != this->_bcValues.end(); ++bc_iter) {
3511: const int& funcNum = bc_iter->first;
3512: const bc_value_type& funcVal = bc_iter->second;
3513: MPI_Bcast((void *) &funcNum, 1, MPI_INT, 0, this->comm());
3514: MPI_Bcast((void *) &funcVal, 4, MPI_DOUBLE, 0, this->comm());
3515: }
3516: }
3517: };
3518: public: // BC values for PyLith
3519: const Obj<foliated_section_type>& getBoundariesNew() {
3520: if (this->_boundaries.isNull()) {
3521: this->_boundaries = new foliated_section_type(this->getTopology());
3522: }
3523: return this->_boundaries;
3524: };
3525: public: // Discretization
3526: void setupField(const Obj<real_section_type>& s, const bool postponeGhosts = false) {
3527: const std::string& name = this->_boundaryCondition->getLabelName();
3528: const patch_type patch = 0;
3530: for(int d = 0; d <= this->_dim; ++d) {
3531: s->setFiberDimensionByDepth(patch, d, this->_discretization->getNumDof(d));
3532: }
3533: if (!name.empty()) {
3534: const Obj<topology_type::label_sequence>& boundary = this->_topology->getLabelStratum(patch, name, 1);
3536: for(topology_type::label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
3537: s->setFiberDimension(patch, *e_iter, -this->_discretization->getNumDof(this->_topology->depth(patch, *e_iter)));
3538: }
3539: }
3540: s->allocate(postponeGhosts);
3541: if (!name.empty()) {
3542: const Obj<real_section_type>& coordinates = this->getRealSection("coordinates");
3543: const Obj<topology_type::label_sequence>& boundary = this->_topology->getLabelStratum(patch, name, 1);
3545: for(topology_type::label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
3546: const real_section_type::value_type *coords = coordinates->restrictPoint(patch, *e_iter);
3547: const PetscScalar value = this->_boundaryCondition->evaluate(coords);
3549: s->updatePointBC(patch, *e_iter, &value);
3550: }
3551: }
3552: };
3553: public:
3554: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
3555: if (comm == MPI_COMM_NULL) {
3556: comm = this->comm();
3557: }
3558: if (name == "") {
3559: PetscPrintf(comm, "viewing a Mesh\n");
3560: } else {
3561: PetscPrintf(comm, "viewing Mesh '%s'\n", name.c_str());
3562: }
3563: this->getTopology()->view("mesh topology", comm);
3564: Obj<std::set<std::string> > sections = this->getRealSections();
3566: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
3567: this->getRealSection(*name)->view(*name);
3568: }
3569: sections = this->getIntSections();
3570: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
3571: this->getIntSection(*name)->view(*name);
3572: }
3573: sections = this->getPairSections();
3574: for(std::set<std::string>::iterator name = sections->begin(); name != sections->end(); ++name) {
3575: this->getPairSection(*name)->view(*name);
3576: }
3577: };
3578: template<typename value_type>
3579: static std::string printMatrix(const std::string& name, const int rows, const int cols, const value_type matrix[], const int rank = -1)
3580: {
3581: ostringstream output;
3582: ostringstream rankStr;
3584: if (rank >= 0) {
3585: rankStr << "[" << rank << "]";
3586: }
3587: output << rankStr.str() << name << " = " << std::endl;
3588: for(int r = 0; r < rows; r++) {
3589: if (r == 0) {
3590: output << rankStr.str() << " /";
3591: } else if (r == rows-1) {
3592: output << rankStr.str() << " \\";
3593: } else {
3594: output << rankStr.str() << " |";
3595: }
3596: for(int c = 0; c < cols; c++) {
3597: output << " " << matrix[r*cols+c];
3598: }
3599: if (r == 0) {
3600: output << " \\" << std::endl;
3601: } else if (r == rows-1) {
3602: output << " /" << std::endl;
3603: } else {
3604: output << " |" << std::endl;
3605: }
3606: }
3607: return output.str();
3608: };
3609: };
3610: } // namespace ALECompat
3611: #endif