Actual source code: CoSieve.hh
1: #ifndef included_ALE_CoSieve_hh
2: #define included_ALE_CoSieve_hh
4: #ifndef included_ALE_Sieve_hh
5: #include <Sieve.hh>
6: #endif
10: #ifdef PETSC_HAVE_CHACO
11: /* Chaco does not have an include file */
14: float *ewgts, float *x, float *y, float *z, char *outassignname,
15: char *outfilename, short *assignment, int architecture, int ndims_tot,
16: int mesh_dims[3], double *goal, int global_method, int local_method,
17: int rqi_flag, int vmax, int ndims, double eigtol, long seed);
20: }
21: #endif
24: namespace ALE {
25: class ParallelObject {
26: protected:
27: int _debug;
28: MPI_Comm _comm;
29: int _commRank;
30: int _commSize;
31: PetscObject _petscObj;
32: void __init(MPI_Comm comm) {
33: static PetscCookie objType = -1;
34: //const char *id_name = ALE::getClassName<T>();
35: const char *id_name = "ParallelObject";
36: PetscErrorCode ierr;
38: if (objType < 0) {
39: PetscLogClassRegister(&objType, id_name);CHKERROR(ierr, "Error in PetscLogClassRegister");
40: }
41: this->_comm = comm;
42: MPI_Comm_rank(this->_comm, &this->_commRank); CHKERROR(ierr, "Error in MPI_Comm_rank");
43: MPI_Comm_size(this->_comm, &this->_commSize); CHKERROR(ierr, "Error in MPI_Comm_size");
44: PetscObjectCreateGeneric(this->_comm, objType, id_name, &this->_petscObj);CHKERROR(ierr, "Error in PetscObjectCreate");
45: //ALE::restoreClassName<T>(id_name);
46: };
47: public:
48: ParallelObject(MPI_Comm comm = PETSC_COMM_SELF, const int debug = 0) : _debug(debug), _petscObj(NULL) {__init(comm);}
49: virtual ~ParallelObject() {
50: if (this->_petscObj) {
51: PetscErrorCode PetscObjectDestroy(this->_petscObj); CHKERROR(ierr, "Failed in PetscObjectDestroy");
52: this->_petscObj = NULL;
53: }
54: };
55: public:
56: int debug() const {return this->_debug;};
57: void setDebug(const int debug) {this->_debug = debug;};
58: MPI_Comm comm() const {return this->_comm;};
59: int commSize() const {return this->_commSize;};
60: int commRank() const {return this->_commRank;}
61: PetscObject petscObj() const {return this->_petscObj;};
62: };
64: namespace New {
65: template<typename Sieve_>
66: class SieveBuilder {
67: public:
68: typedef Sieve_ sieve_type;
69: typedef std::vector<typename sieve_type::point_type> PointArray;
70: public:
71: static void buildHexFaces(Obj<sieve_type> sieve, int dim, std::map<int, int*>& curElement, std::map<int,PointArray>& bdVertices, std::map<int,PointArray>& faces, typename sieve_type::point_type& cell) {
72: int debug = sieve->debug;
74: if (debug > 1) {std::cout << " Building hex faces for boundary of " << cell << " (size " << bdVertices[dim].size() << "), dim " << dim << std::endl;}
75: if (dim > 3) {
76: throw ALE::Exception("Cannot do hexes of dimension greater than three");
77: } else if (dim > 2) {
78: int nodes[24] = {0, 1, 2, 3, 4, 5, 6, 7, 0, 1, 5, 4,
79: 1, 2, 6, 5, 2, 3, 7, 6, 3, 0, 4, 7};
81: for(int b = 0; b < 6; b++) {
82: typename sieve_type::point_type face;
84: for(int c = 0; c < 4; c++) {
85: bdVertices[dim-1].push_back(bdVertices[dim][nodes[b*4+c]]);
86: }
87: if (debug > 1) {std::cout << " boundary hex face " << b << std::endl;}
88: buildHexFaces(sieve, dim-1, curElement, bdVertices, faces, face);
89: if (debug > 1) {std::cout << " added face " << face << std::endl;}
90: faces[dim].push_back(face);
91: }
92: } else if (dim > 1) {
93: int boundarySize = bdVertices[dim].size();
95: for(int b = 0; b < boundarySize; b++) {
96: typename sieve_type::point_type face;
98: for(int c = 0; c < 2; c++) {
99: bdVertices[dim-1].push_back(bdVertices[dim][(b+c)%boundarySize]);
100: }
101: if (debug > 1) {std::cout << " boundary point " << bdVertices[dim][b] << std::endl;}
102: buildHexFaces(sieve, dim-1, curElement, bdVertices, faces, face);
103: if (debug > 1) {std::cout << " added face " << face << std::endl;}
104: faces[dim].push_back(face);
105: }
106: } else {
107: if (debug > 1) {std::cout << " Just set faces to boundary in 1d" << std::endl;}
108: faces[dim].insert(faces[dim].end(), bdVertices[dim].begin(), bdVertices[dim].end());
109: }
110: if (debug > 1) {
111: for(typename PointArray::iterator f_iter = faces[dim].begin(); f_iter != faces[dim].end(); ++f_iter) {
112: std::cout << " face point " << *f_iter << std::endl;
113: }
114: }
115: // We always create the toplevel, so we could short circuit somehow
116: // Should not have to loop here since the meet of just 2 boundary elements is an element
117: typename PointArray::iterator f_itor = faces[dim].begin();
118: const typename sieve_type::point_type& start = *f_itor;
119: const typename sieve_type::point_type& next = *(++f_itor);
120: Obj<typename sieve_type::supportSet> preElement = sieve->nJoin(start, next, 1);
122: if (preElement->size() > 0) {
123: cell = *preElement->begin();
124: if (debug > 1) {std::cout << " Found old cell " << cell << std::endl;}
125: } else {
126: int color = 0;
128: cell = typename sieve_type::point_type((*curElement[dim])++);
129: for(typename PointArray::iterator f_itor = faces[dim].begin(); f_itor != faces[dim].end(); ++f_itor) {
130: sieve->addArrow(*f_itor, cell, color++);
131: }
132: if (debug > 1) {std::cout << " Added cell " << cell << " dim " << dim << std::endl;}
133: }
134: };
135: static void buildFaces(Obj<sieve_type> sieve, int dim, std::map<int, int*>& curElement, std::map<int,PointArray>& bdVertices, std::map<int,PointArray>& faces, typename sieve_type::point_type& cell) {
136: int debug = sieve->debug;
138: if (debug > 1) {
139: if (cell >= 0) {
140: std::cout << " Building faces for boundary of " << cell << " (size " << bdVertices[dim].size() << "), dim " << dim << std::endl;
141: } else {
142: std::cout << " Building faces for boundary of undetermined cell (size " << bdVertices[dim].size() << "), dim " << dim << std::endl;
143: }
144: }
145: faces[dim].clear();
146: if (dim > 1) {
147: // Use the cone construction
148: for(typename PointArray::iterator b_itor = bdVertices[dim].begin(); b_itor != bdVertices[dim].end(); ++b_itor) {
149: typename sieve_type::point_type face = -1;
151: bdVertices[dim-1].clear();
152: for(typename PointArray::iterator i_itor = bdVertices[dim].begin(); i_itor != bdVertices[dim].end(); ++i_itor) {
153: if (i_itor != b_itor) {
154: bdVertices[dim-1].push_back(*i_itor);
155: }
156: }
157: if (debug > 1) {std::cout << " boundary point " << *b_itor << std::endl;}
158: buildFaces(sieve, dim-1, curElement, bdVertices, faces, face);
159: if (debug > 1) {std::cout << " added face " << face << std::endl;}
160: faces[dim].push_back(face);
161: }
162: } else {
163: if (debug > 1) {std::cout << " Just set faces to boundary in 1d" << std::endl;}
164: faces[dim].insert(faces[dim].end(), bdVertices[dim].begin(), bdVertices[dim].end());
165: }
166: if (debug > 1) {
167: for(typename PointArray::iterator f_iter = faces[dim].begin(); f_iter != faces[dim].end(); ++f_iter) {
168: std::cout << " face point " << *f_iter << std::endl;
169: }
170: }
171: // We always create the toplevel, so we could short circuit somehow
172: // Should not have to loop here since the meet of just 2 boundary elements is an element
173: typename PointArray::iterator f_itor = faces[dim].begin();
174: const typename sieve_type::point_type& start = *f_itor;
175: const typename sieve_type::point_type& next = *(++f_itor);
176: Obj<typename sieve_type::supportSet> preElement = sieve->nJoin(start, next, 1);
178: if (preElement->size() > 0) {
179: cell = *preElement->begin();
180: if (debug > 1) {std::cout << " Found old cell " << cell << std::endl;}
181: } else {
182: int color = 0;
184: cell = typename sieve_type::point_type((*curElement[dim])++);
185: for(typename PointArray::iterator f_itor = faces[dim].begin(); f_itor != faces[dim].end(); ++f_itor) {
186: sieve->addArrow(*f_itor, cell, color++);
187: }
188: if (debug > 1) {std::cout << " Added cell " << cell << " dim " << dim << std::endl;}
189: }
190: };
194: // Build a topology from a connectivity description
195: // (0, 0) ... (0, numCells-1): dim-dimensional cells
196: // (0, numCells) ... (0, numVertices): vertices
197: // The other cells are numbered as they are requested
198: static void buildTopology(Obj<sieve_type> sieve, int dim, int numCells, int cells[], int numVertices, bool interpolate = true, int corners = -1) {
199: int debug = sieve->debug;
201: ALE_LOG_EVENT_BEGIN;
202: if (sieve->commRank() != 0) {
203: ALE_LOG_EVENT_END;
204: return;
205: }
206: // Create a map from dimension to the current element number for that dimension
207: std::map<int,int*> curElement;
208: std::map<int,PointArray> bdVertices;
209: std::map<int,PointArray> faces;
210: int curCell = 0;
211: int curVertex = numCells;
212: int newElement = numCells+numVertices;
214: if (corners < 0) corners = dim+1;
215: curElement[0] = &curVertex;
216: curElement[dim] = &curCell;
217: for(int d = 1; d < dim; d++) {
218: curElement[d] = &newElement;
219: }
220: for(int c = 0; c < numCells; c++) {
221: typename sieve_type::point_type cell(c);
223: // Build the cell
224: if (interpolate) {
225: bdVertices[dim].clear();
226: for(int b = 0; b < corners; b++) {
227: typename sieve_type::point_type vertex(cells[c*corners+b]+numCells);
229: if (debug > 1) {std::cout << "Adding boundary vertex " << vertex << std::endl;}
230: bdVertices[dim].push_back(vertex);
231: }
232: if (debug) {std::cout << "cell " << cell << " num boundary vertices " << bdVertices[dim].size() << std::endl;}
234: if (corners != dim+1) {
235: buildHexFaces(sieve, dim, curElement, bdVertices, faces, cell);
236: } else {
237: buildFaces(sieve, dim, curElement, bdVertices, faces, cell);
238: }
239: } else {
240: for(int b = 0; b < corners; b++) {
241: sieve->addArrow(typename sieve_type::point_type(cells[c*corners+b]+numCells), cell, b);
242: }
243: if (debug) {
244: if (debug > 1) {
245: for(int b = 0; b < corners; b++) {
246: std::cout << " Adding vertex " << typename sieve_type::point_type(cells[c*corners+b]+numCells) << std::endl;
247: }
248: }
249: std::cout << "Adding cell " << cell << " dim " << dim << std::endl;
250: }
251: }
252: }
253: ALE_LOG_EVENT_END;
254: };
255: };
257: // A Topology is a collection of Sieves
258: // Each Sieve has a label, which we call a \emph{patch}
259: // The collection itself we call a \emph{sheaf}
260: // The main operation we provide in Topology is the creation of a \emph{label}
261: // A label is a bidirectional mapping of Sieve points to integers, implemented with a Sifter
262: template<typename Patch_, typename Sieve_>
263: class Topology : public ALE::ParallelObject {
264: public:
265: typedef Patch_ patch_type;
266: typedef Sieve_ sieve_type;
267: typedef typename sieve_type::point_type point_type;
268: typedef typename std::map<patch_type, Obj<sieve_type> > sheaf_type;
269: typedef typename ALE::Sifter<int, point_type, int> patch_label_type;
270: typedef typename std::map<patch_type, Obj<patch_label_type> > label_type;
271: typedef typename std::map<patch_type, int> max_label_type;
272: typedef typename std::map<const std::string, label_type> labels_type;
273: typedef typename patch_label_type::supportSequence label_sequence;
274: typedef typename std::set<point_type> point_set_type;
275: protected:
276: sheaf_type _sheaf;
277: labels_type _labels;
278: int _maxHeight;
279: max_label_type _maxHeights;
280: int _maxDepth;
281: max_label_type _maxDepths;
282: // Work space
283: Obj<point_set_type> _modifiedPoints;
284: public:
285: Topology(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug), _maxHeight(-1), _maxDepth(-1) {
286: this->_modifiedPoints = new point_set_type();
287: };
288: public: // Verifiers
289: void checkPatch(const patch_type& patch) {
290: if (this->_sheaf.find(patch) == this->_sheaf.end()) {
291: ostringstream msg;
292: msg << "Invalid topology patch: " << patch << std::endl;
293: throw ALE::Exception(msg.str().c_str());
294: }
295: };
296: void checkLabel(const std::string& name) {
297: if (this->_labels.find(name) == this->_labels.end()) {
298: ostringstream msg;
299: msg << "Invalid label name: " << name << std::endl;
300: throw ALE::Exception(msg.str().c_str());
301: }
302: };
303: public: // Accessors
304: const Obj<sieve_type>& getPatch(const patch_type& patch) {
305: this->checkPatch(patch);
306: return this->_sheaf[patch];
307: };
308: void setPatch(const patch_type& patch, const Obj<sieve_type>& sieve) {
309: this->_sheaf[patch] = sieve;
310: };
311: int getValue (const Obj<patch_label_type>& label, const point_type& point, const int defValue = 0) {
312: const Obj<typename patch_label_type::coneSequence>& cone = label->cone(point);
314: if (cone->size() == 0) return defValue;
315: return *cone->begin();
316: };
317: template<typename InputPoints>
318: int getMaxValue (const Obj<patch_label_type>& label, const Obj<InputPoints>& points, const int defValue = 0) {
319: int maxValue = defValue;
321: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
322: maxValue = std::max(maxValue, this->getValue(label, *p_iter, defValue));
323: }
324: return maxValue;
325: };
326: void setValue(const Obj<patch_label_type>& label, const point_type& point, const int value) {
327: label->setCone(value, point);
328: };
329: const Obj<patch_label_type>& createLabel(const patch_type& patch, const std::string& name) {
330: this->checkPatch(patch);
331: if (this->_labels.find(name) == this->_labels.end()) {
332: this->_labels[name][patch] = new patch_label_type(this->comm(), this->debug());
333: }
334: return this->_labels[name][patch];
335: };
336: const Obj<patch_label_type>& getLabel(const patch_type& patch, const std::string& name) {
337: this->checkLabel(name);
338: this->checkPatch(patch);
339: return this->_labels[name][patch];
340: };
341: const Obj<label_sequence>& getLabelStratum(const patch_type& patch, const std::string& name, int label) {
342: this->checkLabel(name);
343: this->checkPatch(patch);
344: return this->_labels[name][patch]->support(label);
345: };
346: const sheaf_type& getPatches() {
347: return this->_sheaf;
348: };
349: void clear() {
350: this->_sheaf.clear();
351: this->_labels.clear();
352: this->_maxHeight = -1;
353: this->_maxHeights.clear();
354: this->_maxDepth = -1;
355: this->_maxDepths.clear();
356: };
357: public:
358: template<class InputPoints>
359: void computeHeight(const Obj<patch_label_type>& height, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxHeight) {
360: this->_modifiedPoints->clear();
362: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
363: // Compute the max height of the points in the support of p, and add 1
364: int h0 = this->getValue(height, *p_iter, -1);
365: int h1 = this->getMaxValue(height, sieve->support(*p_iter), -1) + 1;
367: if(h1 != h0) {
368: this->setValue(height, *p_iter, h1);
369: if (h1 > maxHeight) maxHeight = h1;
370: this->_modifiedPoints->insert(*p_iter);
371: }
372: }
373: // FIX: We would like to avoid the copy here with cone()
374: if(this->_modifiedPoints->size() > 0) {
375: this->computeHeight(height, sieve, sieve->cone(this->_modifiedPoints), maxHeight);
376: }
377: };
378: void computeHeights() {
379: const std::string name("height");
381: this->_maxHeight = -1;
382: for(typename sheaf_type::iterator s_iter = this->_sheaf.begin(); s_iter != this->_sheaf.end(); ++s_iter) {
383: Obj<patch_label_type> label = this->createLabel(s_iter->first, name);
385: this->_maxHeights[s_iter->first] = -1;
386: this->computeHeight(label, s_iter->second, s_iter->second->leaves(), this->_maxHeights[s_iter->first]);
387: if (this->_maxHeights[s_iter->first] > this->_maxHeight) this->_maxHeight = this->_maxHeights[s_iter->first];
388: }
389: };
390: int height() {return this->_maxHeight;};
391: int height(const patch_type& patch) {
392: this->checkPatch(patch);
393: return this->_maxHeights[patch];
394: };
395: int height(const patch_type& patch, const point_type& point) {
396: return this->getValue(this->_labels["height"][patch], point, -1);
397: };
398: const Obj<label_sequence>& heightStratum(const patch_type& patch, int height) {
399: return this->getLabelStratum(patch, "height", height);
400: };
401: template<class InputPoints>
402: void computeDepth(const Obj<patch_label_type>& depth, const Obj<sieve_type>& sieve, const Obj<InputPoints>& points, int& maxDepth) {
403: this->_modifiedPoints->clear();
405: for(typename InputPoints::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
406: // Compute the max depth of the points in the cone of p, and add 1
407: int d0 = this->getValue(depth, *p_iter, -1);
408: int d1 = this->getMaxValue(depth, sieve->cone(*p_iter), -1) + 1;
410: if(d1 != d0) {
411: this->setValue(depth, *p_iter, d1);
412: if (d1 > maxDepth) maxDepth = d1;
413: this->_modifiedPoints->insert(*p_iter);
414: }
415: }
416: // FIX: We would like to avoid the copy here with support()
417: if(this->_modifiedPoints->size() > 0) {
418: this->computeDepth(depth, sieve, sieve->support(this->_modifiedPoints), maxDepth);
419: }
420: };
421: void computeDepths() {
422: const std::string name("depth");
424: this->_maxDepth = -1;
425: for(typename sheaf_type::iterator s_iter = this->_sheaf.begin(); s_iter != this->_sheaf.end(); ++s_iter) {
426: Obj<patch_label_type> label = this->createLabel(s_iter->first, name);
428: this->_maxDepths[s_iter->first] = -1;
429: this->computeDepth(label, s_iter->second, s_iter->second->roots(), this->_maxDepths[s_iter->first]);
430: if (this->_maxDepths[s_iter->first] > this->_maxDepth) this->_maxDepth = this->_maxDepths[s_iter->first];
431: }
432: };
433: int depth() {return this->_maxDepth;};
434: int depth(const patch_type& patch) {
435: this->checkPatch(patch);
436: return this->_maxDepths[patch];
437: };
438: int depth(const patch_type& patch, const point_type& point) {
439: return this->getValue(this->_labels["depth"][patch], point, -1);
440: };
441: const Obj<label_sequence>& depthStratum(const patch_type& patch, int depth) {
442: return this->getLabelStratum(patch, "depth", depth);
443: };
446: void stratify() {
447: ALE_LOG_EVENT_BEGIN;
448: this->computeHeights();
449: this->computeDepths();
450: ALE_LOG_EVENT_END;
451: };
452: public: // Viewers
453: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
454: if (comm == MPI_COMM_NULL) {
455: comm = this->comm();
456: }
457: if (name == "") {
458: PetscPrintf(comm, "viewing a Topology\n");
459: } else {
460: PetscPrintf(comm, "viewing Topology '%s'\n", name.c_str());
461: }
462: for(typename sheaf_type::const_iterator s_iter = this->_sheaf.begin(); s_iter != this->_sheaf.end(); ++s_iter) {
463: ostringstream txt;
465: txt << "Patch " << s_iter->first;
466: s_iter->second->view(txt.str().c_str(), comm);
467: }
468: };
469: };
471: template<typename Topology_>
472: class Partitioner {
473: public:
474: typedef Topology_ topology_type;
475: typedef typename topology_type::sieve_type sieve_type;
476: typedef typename topology_type::patch_type patch_type;
477: typedef typename topology_type::point_type point_type;
478: public:
481: // This creates a CSR representation of the adjacency matrix for cells
482: static void buildDualCSR(const Obj<topology_type>& topology, const int dim, const patch_type& patch, int **offsets, int **adjacency) {
483: ALE_LOG_EVENT_BEGIN;
484: const Obj<sieve_type>& sieve = topology->getPatch(patch);
485: const Obj<typename topology_type::label_sequence>& elements = topology->heightStratum(patch, 0);
486: int numElements = elements->size();
487: int corners = sieve->cone(*elements->begin())->size();
488: int *off = new int[numElements+1];
490: std::set<point_type> *neighborCells = new std::set<point_type>[numElements];
491: int faceVertices = -1;
493: if (topology->depth(patch) != 1) {
494: throw ALE::Exception("Not yet implemented for interpolated meshes");
495: }
496: if (corners == dim+1) {
497: faceVertices = dim;
498: } else if ((dim == 2) && (corners == 4)) {
499: faceVertices = 2;
500: } else if ((dim == 3) && (corners == 8)) {
501: faceVertices = 4;
502: } else {
503: throw ALE::Exception("Could not determine number of face vertices");
504: }
505: for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
506: const Obj<typename sieve_type::traits::coneSequence>& vertices = sieve->cone(*e_iter);
507: typename sieve_type::traits::coneSequence::iterator vEnd = vertices->end();
509: for(typename sieve_type::traits::coneSequence::iterator v_iter = vertices->begin(); v_iter != vEnd; ++v_iter) {
510: const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*v_iter);
511: typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();
513: for(typename sieve_type::traits::supportSequence::iterator n_iter = neighbors->begin(); n_iter != nEnd; ++n_iter) {
514: if (*e_iter == *n_iter) continue;
515: if ((int) sieve->meet(*e_iter, *n_iter)->size() == faceVertices) {
516: neighborCells[*e_iter].insert(*n_iter);
517: }
518: }
519: }
520: }
521: off[0] = 0;
522: for(int e = 1; e <= numElements; e++) {
523: off[e] = neighborCells[e-1].size() + off[e-1];
524: }
525: int *adj = new int[off[numElements]];
526: int offset = 0;
527: for(int e = 0; e < numElements; e++) {
528: for(typename std::set<point_type>::iterator n_iter = neighborCells[e].begin(); n_iter != neighborCells[e].end(); ++n_iter) {
529: adj[offset++] = *n_iter;
530: }
531: }
532: delete [] neighborCells;
533: if (offset != off[numElements]) {
534: ostringstream msg;
535: msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
536: throw ALE::Exception(msg.str().c_str());
537: }
538: *offsets = off;
539: *adjacency = adj;
540: ALE_LOG_EVENT_END;
541: };
542: #ifdef PETSC_HAVE_CHACO
545: static short *partitionSieve_Chaco(const Obj<topology_type>& topology, const int dim) {
546: short *assignment = NULL; /* set number of each vtx (length n) */
548: ALE_LOG_EVENT_BEGIN;
549: if (topology->commRank() == 0) {
550: /* arguments for Chaco library */
551: FREE_GRAPH = 0; /* Do not let Chaco free my memory */
552: int nvtxs; /* number of vertices in full graph */
553: int *start; /* start of edge list for each vertex */
554: int *adjacency; /* = adj -> j; edge list data */
555: int *vwgts = NULL; /* weights for all vertices */
556: float *ewgts = NULL; /* weights for all edges */
557: float *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
558: char *outassignname = NULL; /* name of assignment output file */
559: char *outfilename = NULL; /* output file name */
560: int architecture = 1; /* 0 => hypercube, d => d-dimensional mesh */
561: int ndims_tot = 0; /* total number of cube dimensions to divide */
562: int mesh_dims[3]; /* dimensions of mesh of processors */
563: double *goal = NULL; /* desired set sizes for each set */
564: int global_method = 1; /* global partitioning algorithm */
565: int local_method = 1; /* local partitioning algorithm */
566: int rqi_flag = 0; /* should I use RQI/Symmlq eigensolver? */
567: int vmax = 200; /* how many vertices to coarsen down to? */
568: int ndims = 1; /* number of eigenvectors (2^d sets) */
569: double eigtol = 0.001; /* tolerance on eigenvectors */
570: long seed = 123636512; /* for random graph mutations */
571: int patch = 0;
574: nvtxs = topology->heightStratum(patch, 0)->size();
575: mesh_dims[0] = topology->commSize(); mesh_dims[1] = 1; mesh_dims[2] = 1;
576: ALE::New::Partitioner<topology_type>::buildDualCSR(topology, dim, patch, &start, &adjacency);
577: for(int e = 0; e < start[nvtxs]; e++) {
578: adjacency[e]++;
579: }
580: assignment = new short int[nvtxs];
581: PetscMemzero(assignment, nvtxs * sizeof(short));CHKERROR(ierr, "Error in PetscMemzero");
583: /* redirect output to buffer: chaco -> msgLog */
584: #ifdef PETSC_HAVE_UNISTD_H
585: char *msgLog;
586: int fd_stdout, fd_pipe[2], count;
588: fd_stdout = dup(1);
589: pipe(fd_pipe);
590: close(1);
591: dup2(fd_pipe[1], 1);
592: msgLog = new char[16284];
593: #endif
595: interface(nvtxs, start, adjacency, vwgts, ewgts, x, y, z,
596: outassignname, outfilename, assignment, architecture, ndims_tot,
597: mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims,
598: eigtol, seed);
600: #ifdef PETSC_HAVE_UNISTD_H
601: int SIZE_LOG = 10000;
603: fflush(stdout);
604: count = read(fd_pipe[0], msgLog, (SIZE_LOG - 1) * sizeof(char));
605: if (count < 0) count = 0;
606: msgLog[count] = 0;
607: close(1);
608: dup2(fd_stdout, 1);
609: close(fd_stdout);
610: close(fd_pipe[0]);
611: close(fd_pipe[1]);
612: if (topology->debug()) {
613: std::cout << msgLog << std::endl;
614: }
615: delete [] msgLog;
616: #endif
617: delete [] adjacency;
618: delete [] start;
619: }
620: ALE_LOG_EVENT_END;
621: return assignment;
622: };
623: #endif
624: };
626: template<typename Topology_, typename Index_>
627: class Atlas : public ALE::ParallelObject {
628: public:
629: typedef Topology_ topology_type;
630: typedef typename topology_type::patch_type patch_type;
631: typedef typename topology_type::sieve_type sieve_type;
632: typedef typename topology_type::point_type point_type;
633: typedef Index_ index_type;
634: typedef std::vector<index_type> IndexArray;
635: typedef typename std::map<point_type, index_type> chart_type;
636: typedef typename std::map<patch_type, chart_type> indices_type;
637: protected:
638: Obj<topology_type> _topology;
639: indices_type _indices;
640: Obj<IndexArray> _array;
641: public:
642: Atlas(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
643: this->_topology = new topology_type(comm, debug);
644: this->_array = new IndexArray();
645: };
646: Atlas(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()), _topology(topology) {
647: this->_array = new IndexArray();
648: };
649: public: // Accessors
650: const Obj<topology_type>& getTopology() const {return this->_topology;};
651: void setTopology(const Obj<topology_type>& topology) {this->_topology = topology;};
652: void copy(const Obj<Atlas>& atlas) {
653: const typename topology_type::sheaf_type& sheaf = atlas->getTopology()->getPatches();
655: for(typename topology_type::sheaf_type::const_iterator s_iter = sheaf.begin(); s_iter != sheaf.end(); ++s_iter) {
656: const chart_type& chart = atlas->getChart(s_iter->first);
658: for(typename chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
659: this->setFiberDimension(s_iter->first, c_iter->first, c_iter->second.index);
660: }
661: }
662: };
663: void copyByDepth(const Obj<Atlas>& atlas) {
664: this->copyByDepth(atlas, atlas->getTopology());
665: };
666: template<typename AtlasType, typename TopologyType>
667: void copyByDepth(const Obj<AtlasType>& atlas, const Obj<TopologyType>& topology) {
668: const typename topology_type::sheaf_type& patches = topology->getPatches();
669: bool *depths = new bool[topology->depth()+1];
671: for(int d = 0; d <= topology->depth(); d++) depths[d] = false;
672: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
673: const patch_type& patch = p_iter->first;
674: const chart_type& chart = atlas->getChart(p_iter->first);
676: for(typename chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
677: const point_type& point = c_iter->first;
678: const int depth = topology->depth(patch, point);
680: if (!depths[depth]) {
681: this->setFiberDimensionByDepth(patch, depth, atlas->getFiberDimension(patch, point));
682: }
683: }
684: }
685: this->orderPatches();
686: }
687: public: // Verifiers
688: void checkPatch(const patch_type& patch) {
689: this->_topology->checkPatch(patch);
690: if (this->_indices.find(patch) == this->_indices.end()) {
691: ostringstream msg;
692: msg << "Invalid atlas patch: " << patch << std::endl;
693: throw ALE::Exception(msg.str().c_str());
694: }
695: };
696: bool hasPatch(const patch_type& patch) {
697: if (this->_indices.find(patch) == this->_indices.end()) return false;
698: return true;
699: }
700: void clear() {
701: this->_indices.clear();
702: };
703: public: // Sizes
704: int const getFiberDimension(const patch_type& patch, const point_type& p) {
705: return this->_indices[patch][p].index;
706: };
707: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
708: this->_indices[patch][p].prefix = -1;
709: this->_indices[patch][p].index = dim;
710: };
711: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
712: if (this->hasPatch(patch) && (this->_indices[patch].find(p) != this->_indices[patch].end())) {
713: this->_indices[patch][p].index += dim;
714: } else {
715: this->setFiberDimension(patch, p, dim);
716: }
717: };
718: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
719: const Obj<typename topology_type::label_sequence>& points = this->_topology->depthStratum(patch, depth);
721: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
722: this->setFiberDimension(patch, *p_iter, dim);
723: }
724: };
725: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
726: const Obj<typename topology_type::label_sequence>& points = this->_topology->heightStratum(patch, height);
728: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
729: this->setFiberDimension(patch, *p_iter, dim);
730: }
731: };
732: void setFiberDimensionByLabel(const patch_type& patch, const std::string& label, int value, int dim) {
733: const Obj<typename topology_type::label_sequence>& points = this->_topology->getLabelStratum(patch, label, value);
735: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
736: this->setFiberDimension(patch, *p_iter, dim);
737: }
738: };
739: int size(const patch_type& patch) {
740: typename chart_type::iterator end = this->_indices[patch].end();
741: int size = 0;
743: for(typename chart_type::iterator c_iter = this->_indices[patch].begin(); c_iter != end; ++c_iter) {
744: size += c_iter->second.index;
745: }
746: return size;
747: };
748: int size(const patch_type& patch, const point_type& p) {
749: this->checkPatch(patch);
750: Obj<typename sieve_type::coneSet> closure = this->_topology->getPatch(patch)->closure(p);
751: typename sieve_type::coneSet::iterator end = closure->end();
752: int size = 0;
754: for(typename sieve_type::coneSet::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
755: size += this->_indices[patch][*c_iter].index;
756: }
757: return size;
758: };
759: void orderPoint(chart_type& chart, const Obj<sieve_type>& sieve, const point_type& point, int& offset) {
760: const Obj<typename sieve_type::coneSequence>& cone = sieve->cone(point);
761: typename sieve_type::coneSequence::iterator end = cone->end();
763: if (chart[point].prefix < 0) {
764: for(typename sieve_type::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
765: if (this->_debug > 1) {std::cout << " Recursing to " << *c_iter << std::endl;}
766: this->orderPoint(chart, sieve, *c_iter, offset);
767: }
768: if (this->_debug > 1) {std::cout << " Ordering point " << point << " at " << offset << std::endl;}
769: chart[point].prefix = offset;
770: offset += chart[point].index;
771: }
772: }
773: void orderPatch(const patch_type& patch, int& offset) {
774: chart_type& chart = this->_indices[patch];
776: for(typename chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
777: if (this->_debug > 1) {std::cout << "Ordering closure of point " << p_iter->first << std::endl;}
778: this->orderPoint(chart, this->_topology->getPatch(patch), p_iter->first, offset);
779: }
780: };
781: void orderPatches() {
782: for(typename indices_type::iterator i_iter = this->_indices.begin(); i_iter != this->_indices.end(); ++i_iter) {
783: if (this->_debug > 1) {std::cout << "Ordering patch " << i_iter->first << std::endl;}
784: int offset = 0;
786: this->orderPatch(i_iter->first, offset);
787: }
788: };
789: void clearIndices() {
790: for(typename indices_type::iterator i_iter = this->_indices.begin(); i_iter != this->_indices.end(); ++i_iter) {
791: this->_indices[i_iter->first].clear();
792: }
793: };
794: public: // Index retrieval
795: const index_type& getIndex(const patch_type& patch, const point_type& p) {
796: this->checkPatch(patch);
797: return this->_indices[patch][p];
798: };
799: template<typename Numbering>
800: const index_type getIndex(const patch_type& patch, const point_type& p, const Obj<Numbering>& numbering) {
801: this->checkPatch(patch);
802: return index_type(numbering->getIndex(p), this->_indices[patch][p].index);
803: };
804: // Want to return a sequence
805: const Obj<IndexArray>& getIndices(const patch_type& patch, const point_type& p, const int level = -1) {
806: this->_array->clear();
808: if (level == 0) {
809: this->_array->push_back(this->getIndex(patch, p));
810: } else if ((level == 1) || (this->_topology->height(patch) == 1)) {
811: const Obj<typename sieve_type::coneSequence>& cone = this->_topology->getPatch(patch)->cone(p);
813: this->_array->push_back(this->getIndex(patch, p));
814: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
815: this->_array->push_back(this->getIndex(patch, *p_iter));
816: }
817: } else if (level == -1) {
818: Obj<typename sieve_type::coneSet> closure = this->_topology->getPatch(patch)->closure(p);
820: for(typename sieve_type::coneSet::iterator p_iter = closure->begin(); p_iter != closure->end(); ++p_iter) {
821: this->_array->push_back(this->getIndex(patch, *p_iter));
822: }
823: } else {
824: Obj<typename sieve_type::coneArray> cone = this->_topology->getPatch(patch)->nCone(p, level);
826: for(typename sieve_type::coneArray::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
827: this->_array->push_back(this->getIndex(patch, *p_iter));
828: }
829: }
830: return this->_array;
831: };
832: template<typename Numbering>
833: const Obj<IndexArray>& getIndices(const patch_type& patch, const point_type& p, const Obj<Numbering>& numbering, const int level = -1) {
834: this->_array->clear();
836: if (level == 0) {
837: this->_array->push_back(this->getIndex(patch, p, numbering));
838: } else if ((level == 1) || (this->_topology->height(patch) == 1)) {
839: const Obj<typename sieve_type::coneSequence>& cone = this->_topology->getPatch(patch)->cone(p);
841: this->_array->push_back(this->getIndex(patch, p, numbering));
842: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
843: this->_array->push_back(this->getIndex(patch, *p_iter, numbering));
844: }
845: } else if (level == -1) {
846: Obj<typename sieve_type::coneSet> closure = this->_topology->getPatch(patch)->closure(p);
848: for(typename sieve_type::coneSet::iterator p_iter = closure->begin(); p_iter != closure->end(); ++p_iter) {
849: this->_array->push_back(this->getIndex(patch, *p_iter, numbering));
850: }
851: } else {
852: Obj<typename sieve_type::coneArray> cone = this->_topology->getPatch(patch)->nCone(p, level);
854: for(typename sieve_type::coneArray::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
855: this->_array->push_back(this->getIndex(patch, *p_iter, numbering));
856: }
857: }
858: return this->_array;
859: };
860: const chart_type& getChart(const patch_type& patch) {
861: return this->_indices[patch];
862: }
863: };
865: // An AbstractSection is a mapping from Sieve points to sets of values
866: // This is our presentation of a section of a fibre bundle,
867: // in which the Topology is the base space, and
868: // the value sets are vectors in the fiber spaces
869: // The main interface to values is through restrict() and update()
870: // This retrieval uses Sieve closure()
871: // We should call these rawRestrict/rawUpdate
872: // The Section must also be able to set/report the dimension of each fiber
873: // for which we use another section we call an \emph{Atlas}
874: // For some storage schemes, we also want offsets to go with these dimensions
875: // We must have some way of limiting the points associated with values
876: // so each section must support a getPatch() call returning the points with associated fibers
877: // I was using the Chart for this
878: // The Section must be able to participate in \emph{completion}
879: // This means restrict to a provided overlap, and exchange in the restricted sections
880: // Completion does not use hierarchy, so we see the Topology as a DiscreteTopology
882: // A ConstantSection is the simplest Section
883: // All fibers are dimension 1
884: // All values are equal to a constant
885: // We need no value storage and no communication for completion
886: template<typename Topology_, typename Value_>
887: class NewConstantSection : public ALE::ParallelObject {
888: public:
889: typedef Topology_ topology_type;
890: typedef typename topology_type::patch_type patch_type;
891: typedef typename topology_type::sieve_type sieve_type;
892: typedef typename topology_type::point_type point_type;
893: typedef std::set<point_type> chart_type;
894: typedef std::map<patch_type, chart_type> atlas_type;
895: typedef Value_ value_type;
896: protected:
897: Obj<topology_type> _topology;
898: atlas_type _atlas;
899: value_type _value;
900: public:
901: NewConstantSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {};
902: NewConstantSection(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()), _topology(topology) {};
903: NewConstantSection(const Obj<topology_type>& topology, const value_type& value) : ParallelObject(topology->comm(), topology->debug()), _topology(topology), _value(value) {};
904: public: // Verifiers
905: void checkPatch(const patch_type& patch) const {
906: this->_topology->checkPatch(patch);
907: if (this->_atlas.find(patch) == this->_atlas.end()) {
908: ostringstream msg;
909: msg << "Invalid atlas patch " << patch << std::endl;
910: throw ALE::Exception(msg.str().c_str());
911: }
912: };
913: void checkDimension(const int& dim) {
914: if (dim != 1) {
915: ostringstream msg;
916: msg << "Invalid fiber dimension " << dim << " must be 1" << std::endl;
917: throw ALE::Exception(msg.str().c_str());
918: }
919: };
920: bool hasPatch(const patch_type& patch) {
921: if (this->_atlas.find(patch) == this->_atlas.end()) {
922: return false;
923: }
924: return true;
925: };
926: public: // Accessors
927: const Obj<topology_type>& getTopology() const {return this->_topology;};
928: void setTopology(const Obj<topology_type>& topology) {this->_topology = topology;};
929: const chart_type& getPatch(const patch_type& patch) {
930: this->checkPatch(patch);
931: return this->_atlas[patch];
932: };
933: void updatePatch(const patch_type& patch, const point_type& point) {
934: this->_atlas[patch].insert(point);
935: };
936: template<typename Points>
937: void updatePatch(const patch_type& patch, const Obj<Points>& points) {
938: this->_atlas[patch].insert(points->begin(), points->end());
939: };
940: public: // Sizes
941: void clear() {
942: this->_atlas.clear();
943: };
944: int getFiberDimension(const patch_type& patch, const point_type& p) const {return 1;};
945: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
946: this->checkDimension(dim);
947: this->updatePatch(patch, p);
948: };
949: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
950: if (this->hasPatch(patch) && (this->_atlas[patch].find(p) != this->_atlas[patch].end())) {
951: ostringstream msg;
952: msg << "Invalid addition to fiber dimension " << dim << " cannot exceed 1" << std::endl;
953: throw ALE::Exception(msg.str().c_str());
954: } else {
955: this->setFiberDimension(patch, p, dim);
956: }
957: };
958: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
959: this->setFiberDimensionByLabel(patch, "depth", depth, dim);
960: };
961: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
962: this->setFiberDimensionByLabel(patch, "height", height, dim);
963: };
964: void setFiberDimensionByLabel(const patch_type& patch, const std::string& label, int value, int dim) {
965: const Obj<typename topology_type::label_sequence>& points = this->_topology->getLabelStratum(patch, label, value);
967: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
968: this->setFiberDimension(patch, *p_iter, dim);
969: }
970: };
971: int size(const patch_type& patch) {return this->_atlas[patch].size();};
972: int size(const patch_type& patch, const point_type& p) {return 1;};
973: void orderPatches() {};
974: public: // Restriction
975: const value_type *restrict(const patch_type& patch, const point_type& p) const {
976: this->checkPatch(patch);
977: return &this->_value;
978: };
979: const value_type *restrictPoint(const patch_type& patch, const point_type& p) const {return this->restrict(patch, p);};
980: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
981: this->checkPatch(patch);
982: this->_value = v[0];
983: };
984: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {return this->update(patch, p, v);};
985: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
986: this->checkPatch(patch);
987: this->_value += v[0];
988: };
989: public:
990: void copy(const NewConstantSection& section) {
991: this->_value = section->restrict(this->_topology->getPatches().begin().first, point_type());
992: };
993: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
994: ostringstream txt;
995: int rank;
997: if (comm == MPI_COMM_NULL) {
998: comm = this->comm();
999: rank = this->commRank();
1000: } else {
1001: MPI_Comm_rank(comm, &rank);
1002: }
1003: if (name == "") {
1004: if(rank == 0) {
1005: txt << "viewing a ConstantSection" << std::endl;
1006: }
1007: } else {
1008: if(rank == 0) {
1009: txt << "viewing ConstantSection '" << name << "'" << std::endl;
1010: }
1011: }
1012: const typename topology_type::sheaf_type& sheaf = this->_topology->getPatches();
1014: for(typename topology_type::sheaf_type::const_iterator p_iter = sheaf.begin(); p_iter != sheaf.end(); ++p_iter) {
1015: txt <<"["<<this->commRank()<<"]: Patch " << p_iter->first << std::endl;
1016: txt <<"["<<this->commRank()<<"]: Value " << this->_value << std::endl;
1017: }
1018: PetscSynchronizedPrintf(comm, txt.str().c_str());
1019: PetscSynchronizedFlush(comm);
1020: };
1021: };
1023: // A UniformSection often acts as an Atlas
1024: // All fibers are the same dimension
1025: // Note we can use a ConstantSection for this Atlas
1026: // Each point may have a different vector
1027: // Thus we need storage for values, and hence must implement completion
1028: template<typename Topology_, typename Value_, int fiberDim = 1>
1029: class UniformSection : public ALE::ParallelObject {
1030: public:
1031: typedef Topology_ topology_type;
1032: typedef typename topology_type::patch_type patch_type;
1033: typedef typename topology_type::sieve_type sieve_type;
1034: typedef typename topology_type::point_type point_type;
1035: typedef NewConstantSection<topology_type, int> atlas_type;
1036: typedef typename atlas_type::chart_type chart_type;
1037: typedef Value_ value_type;
1038: typedef struct {value_type v[fiberDim];} fiber_type;
1039: typedef std::map<point_type, fiber_type> array_type;
1040: typedef std::map<patch_type, array_type> values_type;
1041: protected:
1042: Obj<atlas_type> _atlas;
1043: values_type _arrays;
1044: public:
1045: UniformSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
1046: this->_atlas = new atlas_type(comm, debug);
1047: };
1048: UniformSection(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()) {
1049: this->_atlas = new atlas_type(topology);
1050: };
1051: UniformSection(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas) {};
1052: protected:
1053: value_type *getRawArray(const int size) {
1054: static value_type *array = NULL;
1055: static int maxSize = 0;
1057: if (size > maxSize) {
1058: maxSize = size;
1059: if (array) delete [] array;
1060: array = new value_type[maxSize];
1061: };
1062: return array;
1063: };
1064: public: // Verifiers
1065: void checkPatch(const patch_type& patch) {
1066: this->_atlas->checkPatch(patch);
1067: if (this->_arrays.find(patch) == this->_arrays.end()) {
1068: ostringstream msg;
1069: msg << "Invalid section patch: " << patch << std::endl;
1070: throw ALE::Exception(msg.str().c_str());
1071: }
1072: };
1073: void checkDimension(const int& dim) {
1074: if (dim != fiberDim) {
1075: ostringstream msg;
1076: msg << "Invalid fiber dimension " << dim << " must be " << fiberDim << std::endl;
1077: throw ALE::Exception(msg.str().c_str());
1078: }
1079: };
1080: public: // Accessors
1081: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
1082: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
1083: const Obj<topology_type>& getTopology() {return this->_atlas->getTopology();};
1084: void setTopology(const Obj<topology_type>& topology) {this->_atlas->setTopology(topology);};
1085: const chart_type& getPatch(const patch_type& patch) {
1086: return this->_atlas->getPatch(patch);
1087: };
1088: void updatePatch(const patch_type& patch, const point_type& point) {
1089: this->setFiberDimension(patch, point, 1);
1090: };
1091: template<typename Points>
1092: void updatePatch(const patch_type& patch, const Obj<Points>& points) {
1093: for(typename Points::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
1094: this->setFiberDimension(patch, *p_iter, 1);
1095: }
1096: };
1097: public: // Sizes
1098: void clear() {
1099: this->_atlas->clear();
1100: this->_arrays.clear();
1101: };
1102: int getFiberDimension(const patch_type& patch, const point_type& p) const {
1103: // Could check for non-existence here
1104: return this->_atlas->restrict(patch, p)[0];
1105: };
1106: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1107: this->checkDimension(dim);
1108: this->_atlas->updatePatch(patch, p);
1109: this->_atlas->update(patch, p, &dim);
1110: };
1111: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1112: if (this->hasPatch(patch) && (this->_atlas[patch].find(p) != this->_atlas[patch].end())) {
1113: ostringstream msg;
1114: msg << "Invalid addition to fiber dimension " << dim << " cannot exceed " << fiberDim << std::endl;
1115: throw ALE::Exception(msg.str().c_str());
1116: } else {
1117: this->setFiberDimension(patch, p, dim);
1118: }
1119: };
1120: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
1121: this->setFiberDimensionByLabel(patch, "depth", depth, dim);
1122: };
1123: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
1124: this->setFiberDimensionByLabel(patch, "height", height, dim);
1125: };
1126: void setFiberDimensionByLabel(const patch_type& patch, const std::string& label, int value, int dim) {
1127: const Obj<typename topology_type::label_sequence>& points = this->getTopology()->getLabelStratum(patch, label, value);
1129: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
1130: this->setFiberDimension(patch, *p_iter, dim);
1131: }
1132: };
1133: int size(const patch_type& patch) {
1134: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
1135: int size = 0;
1137: for(typename atlas_type::chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
1138: size += this->getFiberDimension(patch, *p_iter);
1139: }
1140: return size;
1141: };
1142: int size(const patch_type& patch, const point_type& p) {
1143: Obj<typename sieve_type::coneSet> closure = this->getTopology()->getPatch(patch)->closure(p);
1144: typename sieve_type::coneSet::iterator end = closure->end();
1145: int size = 0;
1147: for(typename sieve_type::coneSet::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
1148: size += this->getFiberDimension(patch, *c_iter);
1149: }
1150: return size;
1151: };
1152: void orderPatches() {};
1153: public: // Restriction
1154: const array_type& restrict(const patch_type& patch) {
1155: this->checkPatch(patch);
1156: return this->_arrays[patch];
1157: };
1158: // Return a smart pointer?
1159: const value_type *restrict(const patch_type& patch, const point_type& p) {
1160: this->checkPatch(patch);
1161: array_type& array = this->_arrays[patch];
1162: value_type *values = this->getRawArray(this->size(patch, p));
1164: if (this->getTopology()->height(patch) == 1) {
1165: // Only avoids the copy of closure()
1166: const int& dim = this->_atlas->restrict(patch, p)[0];
1167: int j = -1;
1169: for(int i = 0; i < dim; ++i) {
1170: values[++j] = array[p].v[i];
1171: }
1172: // Should be closure()
1173: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1175: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1176: const int& dim = this->_atlas->restrict(patch, *p_iter)[0];
1178: for(int i = 0; i < dim; ++i) {
1179: values[++j] = array[*p_iter].v[i];
1180: }
1181: }
1182: } else {
1183: #if 0
1184: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1185: int j = -1;
1187: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1188: const int start = i_iter->prefix;
1189: const int length = i_iter->index;
1191: for(int i = start; i < start + length; ++i) {
1192: values[++j] = a[i];
1193: }
1194: }
1195: #else
1196: throw ALE::Exception("Not yet implemented for interpolated sieves");
1197: #endif
1198: }
1199: return values;
1200: };
1201: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {
1202: this->checkPatch(patch);
1203: return this->_arrays[patch][p].v;
1204: };
1205: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
1206: this->_atlas->checkPatch(patch);
1207: array_type& array = this->_arrays[patch];
1209: if (this->getTopology()->height(patch) == 1) {
1210: // Only avoids the copy of closure()
1211: const int& dim = this->_atlas->restrict(patch, p)[0];
1212: int j = -1;
1214: for(int i = 0; i < dim; ++i) {
1215: array[p].v[i] = v[++j];
1216: }
1217: // Should be closure()
1218: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1220: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1221: const int& dim = this->_atlas->restrict(patch, *p_iter)[0];
1223: for(int i = 0; i < dim; ++i) {
1224: array[*p_iter].v[i] = v[++j];
1225: }
1226: }
1227: } else {
1228: #if 0
1229: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1230: int j = -1;
1232: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1233: const int start = i_iter->prefix;
1234: const int length = i_iter->index;
1236: for(int i = start; i < start + length; ++i) {
1237: a[i] = v[++j];
1238: }
1239: }
1240: #else
1241: throw ALE::Exception("Not yet implemented for interpolated sieves");
1242: #endif
1243: }
1244: };
1245: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
1246: this->_atlas->checkPatch(patch);
1247: array_type& array = this->_arrays[patch];
1249: if (this->getTopology()->height(patch) == 1) {
1250: // Only avoids the copy of closure()
1251: const int& dim = this->_atlas->restrict(patch, p)[0];
1252: int j = -1;
1254: for(int i = 0; i < dim; ++i) {
1255: array[p].v[i] += v[++j];
1256: }
1257: // Should be closure()
1258: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1260: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1261: const int& dim = this->_atlas->restrict(patch, *p_iter)[0];
1263: for(int i = 0; i < dim; ++i) {
1264: array[*p_iter].v[i] += v[++j];
1265: }
1266: }
1267: } else {
1268: throw ALE::Exception("Not yet implemented for interpolated sieves");
1269: }
1270: };
1271: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
1272: this->_atlas->checkPatch(patch);
1273: for(int i = 0; i < fiberDim; ++i) {
1274: this->_arrays[patch][p].v[i] = v[i];
1275: }
1276: };
1277: public:
1278: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) {
1279: ostringstream txt;
1280: int rank;
1282: if (comm == MPI_COMM_NULL) {
1283: comm = this->comm();
1284: rank = this->commRank();
1285: } else {
1286: MPI_Comm_rank(comm, &rank);
1287: }
1288: if (name == "") {
1289: if(rank == 0) {
1290: txt << "viewing a UniformSection" << std::endl;
1291: }
1292: } else {
1293: if(rank == 0) {
1294: txt << "viewing UniformSection '" << name << "'" << std::endl;
1295: }
1296: }
1297: for(typename values_type::const_iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
1298: const patch_type& patch = a_iter->first;
1299: array_type& array = this->_arrays[patch];
1301: txt << "[" << this->commRank() << "]: Patch " << patch << std::endl;
1302: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
1304: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
1305: const point_type& point = *p_iter;
1306: const typename atlas_type::value_type dim = this->_atlas->restrict(patch, point)[0];
1308: if (dim != 0) {
1309: txt << "[" << this->commRank() << "]: " << point << " dim " << dim << " ";
1310: for(int i = 0; i < dim; i++) {
1311: txt << " " << array[point].v[i];
1312: }
1313: txt << std::endl;
1314: }
1315: }
1316: }
1317: PetscSynchronizedPrintf(comm, txt.str().c_str());
1318: PetscSynchronizedFlush(comm);
1319: };
1320: };
1322: // A Numbering is a UniformSection giving each Sieve point a consecutive number
1323: // We need to support versions which operate both intra- and inter-patch
1324: // We sometimes need an inverse mapping, so those should really be Sifters or labels
1325: // However, I have not figured out/written distribution for Sifters (do soon)
1327: // A GlobalOrder is an Atlas which calculates interpatch offsets
1329: // A Section is our most general construct (more general ones could be envisioned)
1330: // The Atlas is a UniformSection of dimension 1 and value type Point
1331: // to hold each fiber dimension and offsets into a contiguous patch array
1332: template<typename Topology_, typename Value_>
1333: class NewSection : public ALE::ParallelObject {
1334: public:
1335: typedef Topology_ topology_type;
1336: typedef typename topology_type::patch_type patch_type;
1337: typedef typename topology_type::sieve_type sieve_type;
1338: typedef typename topology_type::point_type point_type;
1339: typedef ALE::Point index_type;
1340: typedef UniformSection<topology_type, index_type> atlas_type;
1341: typedef Value_ value_type;
1342: typedef value_type * array_type;
1343: typedef std::map<patch_type, array_type> values_type;
1344: protected:
1345: Obj<atlas_type> _atlas;
1346: Obj<atlas_type> _atlasNew;
1347: values_type _arrays;
1348: public:
1349: NewSection(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
1350: this->_atlas = new atlas_type(comm, debug);
1351: this->_atlasNew = NULL;
1352: };
1353: NewSection(const Obj<topology_type>& topology) : ParallelObject(topology->comm(), topology->debug()), _atlasNew(NULL) {
1354: this->_atlas = new atlas_type(topology);
1355: };
1356: NewSection(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas), _atlasNew(NULL) {};
1357: virtual ~NewSection() {
1358: for(typename values_type::iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
1359: delete [] a_iter->second;
1360: a_iter->second = NULL;
1361: }
1362: };
1363: protected:
1364: value_type *getRawArray(const int size) {
1365: static value_type *array = NULL;
1366: static int maxSize = 0;
1368: if (size > maxSize) {
1369: maxSize = size;
1370: if (array) delete [] array;
1371: array = new value_type[maxSize];
1372: };
1373: return array;
1374: };
1375: public: // Verifiers
1376: void checkPatch(const patch_type& patch) {
1377: this->_atlas->checkPatch(patch);
1378: if (this->_arrays.find(patch) == this->_arrays.end()) {
1379: ostringstream msg;
1380: msg << "Invalid section patch: " << patch << std::endl;
1381: throw ALE::Exception(msg.str().c_str());
1382: }
1383: };
1384: public: // Accessors
1385: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
1386: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
1387: const Obj<topology_type>& getTopology() {return this->_atlas->getTopology();};
1388: void setTopology(const Obj<topology_type>& topology) {this->_atlas->setTopology(topology);};
1389: public: // Sizes
1390: void clear() {
1391: this->_atlas->clear();
1392: this->_arrays.clear();
1393: };
1394: int getFiberDimension(const patch_type& patch, const point_type& p) const {
1395: // Could check for non-existence here
1396: return this->_atlas->restrict(patch, p)->prefix;
1397: };
1398: void setFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1399: const index_type idx(dim, -1);
1400: this->_atlas->updatePatch(patch, p);
1401: this->_atlas->update(patch, p, &idx);
1402: };
1403: void addFiberDimension(const patch_type& patch, const point_type& p, int dim) {
1404: const value_type values[2] = {dim, 0};
1405: this->_atlas->updatePatch(patch, p);
1406: this->_atlas->updateAdd(patch, p, values);
1407: };
1408: void setFiberDimensionByDepth(const patch_type& patch, int depth, int dim) {
1409: this->setFiberDimensionByLabel(patch, "depth", depth, dim);
1410: };
1411: void setFiberDimensionByHeight(const patch_type& patch, int height, int dim) {
1412: this->setFiberDimensionByLabel(patch, "height", height, dim);
1413: };
1414: void setFiberDimensionByLabel(const patch_type& patch, const std::string& label, int value, int dim) {
1415: const Obj<typename topology_type::label_sequence>& points = this->getTopology()->getLabelStratum(patch, label, value);
1417: for(typename topology_type::label_sequence::iterator p_iter = points->begin(); p_iter != points->end(); ++p_iter) {
1418: this->setFiberDimension(patch, *p_iter, dim);
1419: }
1420: };
1421: int size(const patch_type& patch) {
1422: const typename atlas_type::chart_type& points = this->_atlas->getPatch(patch);
1423: int size = 0;
1425: for(typename atlas_type::chart_type::iterator p_iter = points.begin(); p_iter != points.end(); ++p_iter) {
1426: size += this->getFiberDimension(patch, *p_iter);
1427: }
1428: return size;
1429: };
1430: int size(const patch_type& patch, const point_type& p) {
1431: Obj<typename sieve_type::coneSet> closure = this->getTopology()->getPatch(patch)->closure(p);
1432: typename sieve_type::coneSet::iterator end = closure->end();
1433: int size = 0;
1435: for(typename sieve_type::coneSet::iterator c_iter = closure->begin(); c_iter != end; ++c_iter) {
1436: size += this->getFiberDimension(patch, *c_iter);
1437: }
1438: return size;
1439: };
1440: public: // Allocation
1441: void orderPoint(const Obj<atlas_type>& atlas, const Obj<sieve_type>& sieve, const patch_type& patch, const point_type& point, int& offset) {
1442: const Obj<typename sieve_type::coneSequence>& cone = sieve->cone(point);
1443: typename sieve_type::coneSequence::iterator end = cone->end();
1444: const index_type& idx = atlas->restrict(patch, point)[0];
1445: const int dim = idx.prefix;
1446: index_type newIdx(dim, offset);
1448: if (idx.index < 0) {
1449: for(typename sieve_type::coneSequence::iterator c_iter = cone->begin(); c_iter != end; ++c_iter) {
1450: if (this->_debug > 1) {std::cout << " Recursing to " << *c_iter << std::endl;}
1451: this->orderPoint(atlas, sieve, patch, *c_iter, offset);
1452: }
1453: if (this->_debug > 1) {std::cout << " Ordering point " << point << " at " << offset << std::endl;}
1454: atlas->update(patch, point, &newIdx);
1455: offset += dim;
1456: }
1457: }
1458: void orderPatch(const patch_type& patch, int& offset) {
1459: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
1461: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
1462: if (this->_debug > 1) {std::cout << "Ordering closure of point " << *p_iter << std::endl;}
1463: this->orderPoint(this->_atlas, this->getTopology()->getPatch(patch), patch, *p_iter, offset);
1464: }
1465: };
1466: void orderPatches() {
1467: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
1469: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1470: if (this->_debug > 1) {std::cout << "Ordering patch " << p_iter->first << std::endl;}
1471: int offset = 0;
1473: this->orderPatch(p_iter->first, offset);
1474: }
1475: };
1476: void allocate() {
1477: this->orderPatches();
1478: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
1480: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1481: this->_arrays[p_iter->first] = new value_type[this->size(p_iter->first)];
1482: PetscMemzero(this->_arrays[p_iter->first], this->size(p_iter->first) * sizeof(value_type));
1483: }
1484: };
1485: void addPoint(const patch_type& patch, const point_type& point, const int size) {
1486: const typename atlas_type::chart_type& chart = this->_atlas->getChart(patch);
1488: if (chart.find(point) == chart.end()) {
1489: if (this->_atlasNew.isNull()) {
1490: this->_atlasNew = new atlas_type(this->getTopology());
1491: this->_atlasNew->copy(this->_atlas);
1492: }
1493: this->_atlasNew->setFiberDimension(patch, point, size);
1494: }
1495: };
1496: void reallocate() {
1497: if (this->_atlasNew.isNull()) return;
1498: this->_atlasNew->orderPatches();
1499: const typename topology_type::sheaf_type& patches = this->getTopology()->getPatches();
1501: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1502: const patch_type& patch = p_iter->first;
1503: const typename atlas_type::chart_type& chart = this->_atlas->getChart(patch);
1504: const value_type *array = this->_arrays[patch];
1505: value_type *newArray = new value_type[this->size(patch)];
1507: for(typename atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1508: const int& size = c_iter->second.index;
1509: const int& offset = c_iter->second.prefix;
1510: const int& newOffset = this->_atlasNew->getIndex(patch, c_iter->first).prefix;
1512: for(int i = 0; i < size; ++i) {
1513: newArray[newOffset+i] = array[offset+i];
1514: }
1515: }
1516: delete [] this->_arrays[patch];
1517: this->_arrays[patch] = newArray;
1518: }
1519: this->_atlas = this->_atlasNew;
1520: this->_atlasNew = NULL;
1521: };
1522: public: // Restriction
1523: const value_type *restrict(const patch_type& patch) {
1524: this->checkPatch(patch);
1525: return this->_arrays[patch];
1526: };
1527: // Return a smart pointer?
1528: const value_type *restrict(const patch_type& patch, const point_type& p) {
1529: this->checkPatch(patch);
1530: const value_type *a = this->_arrays[patch];
1531: value_type *values = this->getRawArray(this->size(patch, p));
1533: if (this->getTopology()->height(patch) == 1) {
1534: // Only avoids the copy
1535: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1536: typename sieve_type::coneSequence::iterator end = cone->end();
1537: const index_type& pInd = this->_atlas->restrict(patch, p)[0];
1538: int j = -1;
1540: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
1541: values[++j] = a[i];
1542: }
1543: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
1544: const index_type& ind = this->_atlas->restrict(patch, *p_iter)[0];
1545: const int& start = ind.index;
1546: const int& length = ind.prefix;
1548: for(int i = start; i < start + length; ++i) {
1549: values[++j] = a[i];
1550: }
1551: }
1552: } else {
1553: #if 0
1554: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1555: int j = -1;
1557: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1558: const int start = i_iter->index;
1559: const int length = i_iter->prefix;
1561: for(int i = start; i < start + length; ++i) {
1562: values[++j] = a[i];
1563: }
1564: }
1565: #else
1566: throw ALE::Exception("Not yet implemented for interpolated sieves");
1567: #endif
1568: }
1569: return values;
1570: };
1571: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {
1572: this->checkPatch(patch);
1573: // Using the index structure explicitly
1574: return &(this->_arrays[patch][this->_atlas->restrict(patch, p).index]);
1575: };
1576: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
1577: this->checkPatch(patch);
1578: value_type *a = this->_arrays[patch];
1580: if (this->getTopology()->height(patch) == 1) {
1581: // Only avoids the copy
1582: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1583: typename sieve_type::coneSequence::iterator end = cone->end();
1584: const index_type& pInd = this->_atlas->restrict(patch, p);
1585: int j = -1;
1587: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
1588: a[i] = v[++j];
1589: }
1590: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
1591: const index_type& ind = this->_atlas->restrict(patch, *p_iter);
1592: const int& start = ind.index;
1593: const int& length = ind.prefix;
1595: for(int i = start; i < start + length; ++i) {
1596: a[i] = v[++j];
1597: }
1598: }
1599: } else {
1600: #if 0
1601: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1602: int j = -1;
1604: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1605: const int start = i_iter->index;
1606: const int length = i_iter->prefix;
1608: for(int i = start; i < start + length; ++i) {
1609: a[i] = v[++j];
1610: }
1611: }
1612: #else
1613: throw ALE::Exception("Not yet implemented for interpolated sieves");
1614: #endif
1615: }
1616: };
1617: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
1618: this->checkPatch(patch);
1619: value_type *a = this->_arrays[patch];
1621: if (this->getTopology()->height(patch) == 1) {
1622: // Only avoids the copy
1623: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1624: typename sieve_type::coneSequence::iterator end = cone->end();
1625: const index_type& pInd = this->_atlas->restrict(patch, p)[0];
1626: int j = -1;
1628: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
1629: a[i] += v[++j];
1630: }
1631: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
1632: const index_type& ind = this->_atlas->restrict(patch, *p_iter)[0];
1633: const int& start = ind.index;
1634: const int& length = ind.prefix;
1636: for(int i = start; i < start + length; ++i) {
1637: a[i] += v[++j];
1638: }
1639: }
1640: } else {
1641: #if 0
1642: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1643: int j = -1;
1645: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1646: const int start = i_iter->index;
1647: const int length = i_iter->prefix;
1649: for(int i = start; i < start + length; ++i) {
1650: a[i] += v[++j];
1651: }
1652: }
1653: #else
1654: throw ALE::Exception("Not yet implemented for interpolated sieves");
1655: #endif
1656: }
1657: };
1658: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
1659: this->checkPatch(patch);
1660: const index_type& ind = this->_atlas->getIndex(patch, p);
1661: value_type *a = &(this->_arrays[patch][ind.prefix]);
1663: // Using the index structure explicitly
1664: for(int i = 0; i < ind.index; ++i) {
1665: a[i] = v[i];
1666: }
1667: };
1668: template<typename Input>
1669: void update(const patch_type& patch, const point_type& p, const Obj<Input>& v) {
1670: this->checkPatch(patch);
1671: value_type *a = this->_arrays[patch];
1673: if (this->getTopology()->height(patch) == 1) {
1674: // Only avoids the copy
1675: const Obj<typename sieve_type::coneSequence>& cone = this->getTopology()->getPatch(patch)->cone(p);
1676: typename sieve_type::coneSequence::iterator end = cone->end();
1677: const index_type& pInd = this->_atlas->getIndex(patch, p);
1678: typename Input::iterator v_iter = v->begin();
1679: typename Input::iterator v_end = v->end();
1681: for(int i = pInd.index; i < pInd.prefix + pInd.index; ++i) {
1682: a[i] = *v_iter;
1683: ++v_iter;
1684: }
1685: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != end; ++p_iter) {
1686: const index_type& ind = this->_atlas->getIndex(patch, *p_iter);
1687: const int& start = ind.index;
1688: const int& length = ind.prefix;
1690: for(int i = start; i < start + length; ++i) {
1691: a[i] = *v_iter;
1692: ++v_iter;
1693: }
1694: }
1695: } else {
1696: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1697: typename Input::iterator v_iter = v->begin();
1698: typename Input::iterator v_end = v->end();
1700: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1701: const int& start = i_iter->index;
1702: const int& length = i_iter->prefix;
1704: for(int i = start; i < start + length; ++i) {
1705: a[i] = *v_iter;
1706: ++v_iter;
1707: }
1708: }
1709: }
1710: };
1711: public:
1712: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
1713: ostringstream txt;
1714: int rank;
1716: if (comm == MPI_COMM_NULL) {
1717: comm = this->comm();
1718: rank = this->commRank();
1719: } else {
1720: MPI_Comm_rank(comm, &rank);
1721: }
1722: if (name == "") {
1723: if(rank == 0) {
1724: txt << "viewing a Section" << std::endl;
1725: }
1726: } else {
1727: if(rank == 0) {
1728: txt << "viewing Section '" << name << "'" << std::endl;
1729: }
1730: }
1731: for(typename values_type::const_iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
1732: const patch_type& patch = a_iter->first;
1733: const value_type *array = a_iter->second;
1735: txt << "[" << this->commRank() << "]: Patch " << patch << std::endl;
1736: const typename atlas_type::chart_type& chart = this->_atlas->getPatch(patch);
1738: for(typename atlas_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
1739: const point_type& point = *p_iter;
1740: const index_type& idx = this->_atlas->restrict(patch, point)[0];
1742: if (idx.prefix != 0) {
1743: txt << "[" << this->commRank() << "]: " << point << " dim " << idx.prefix << " offset " << idx.index << " ";
1744: for(int i = 0; i < idx.prefix; i++) {
1745: txt << " " << array[idx.index+i];
1746: }
1747: txt << std::endl;
1748: }
1749: }
1750: }
1751: PetscSynchronizedPrintf(comm, txt.str().c_str());
1752: PetscSynchronizedFlush(comm);
1753: };
1754: };
1756: // An Overlap is a Sifter describing the overlap of two Sieves
1757: // Each arrow is local point ---(remote point)---> remote rank right now
1758: // For XSifter, this should change to (local patch, local point) ---> (remote rank, remote patch, remote point)
1760: template<typename Atlas_, typename Value_>
1761: class Section : public ALE::ParallelObject {
1762: public:
1763: typedef Atlas_ atlas_type;
1764: typedef typename atlas_type::patch_type patch_type;
1765: typedef typename atlas_type::sieve_type sieve_type;
1766: typedef typename atlas_type::topology_type topology_type;
1767: typedef typename atlas_type::point_type point_type;
1768: typedef typename atlas_type::index_type index_type;
1769: typedef Value_ value_type;
1770: typedef std::map<patch_type, value_type *> values_type;
1771: protected:
1772: Obj<atlas_type> _atlas;
1773: Obj<atlas_type> _atlasNew;
1774: values_type _arrays;
1775: public:
1776: Section(MPI_Comm comm, const int debug = 0) : ParallelObject(comm, debug) {
1777: this->_atlas = new atlas_type(comm, debug);
1778: this->_atlasNew = NULL;
1779: };
1780: Section(const Obj<atlas_type>& atlas) : ParallelObject(atlas->comm(), atlas->debug()), _atlas(atlas), _atlasNew(NULL) {};
1781: virtual ~Section() {
1782: for(typename values_type::iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
1783: delete [] a_iter->second;
1784: a_iter->second = NULL;
1785: }
1786: };
1787: protected:
1788: value_type *getRawArray(const int size) {
1789: static value_type *array = NULL;
1790: static int maxSize = 0;
1792: if (size > maxSize) {
1793: maxSize = size;
1794: if (array) delete [] array;
1795: array = new value_type[maxSize];
1796: };
1797: return array;
1798: };
1799: public: // Verifiers
1800: void checkPatch(const patch_type& patch) {
1801: this->_atlas->checkPatch(patch);
1802: if (this->_arrays.find(patch) == this->_arrays.end()) {
1803: ostringstream msg;
1804: msg << "Invalid section patch: " << patch << std::endl;
1805: throw ALE::Exception(msg.str().c_str());
1806: }
1807: };
1808: public: // Accessors
1809: const Obj<atlas_type>& getAtlas() {return this->_atlas;};
1810: void setAtlas(const Obj<atlas_type>& atlas) {this->_atlas = atlas;};
1811: public: // Allocation
1812: void allocate() {
1813: const typename atlas_type::topology_type::sheaf_type& patches = this->_atlas->getTopology()->getPatches();
1815: for(typename atlas_type::topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
1816: this->_arrays[p_iter->first] = new value_type[this->_atlas->size(p_iter->first)];
1817: PetscMemzero(this->_arrays[p_iter->first], this->_atlas->size(p_iter->first) * sizeof(value_type));
1818: }
1819: };
1820: public: // Restriction
1821: const value_type *restrict(const patch_type& patch) {
1822: this->checkPatch(patch);
1823: return this->_arrays[patch];
1824: };
1825: // Return a smart pointer?
1826: const value_type *restrict(const patch_type& patch, const point_type& p) {
1827: this->checkPatch(patch);
1828: const value_type *a = this->_arrays[patch];
1829: value_type *values = this->getRawArray(this->_atlas->size(patch, p));
1831: if (this->_atlas->getTopology()->height(patch) == 1) {
1832: // Only avoids the copy
1833: const Obj<typename sieve_type::coneSequence>& cone = this->_atlas->getTopology()->getPatch(patch)->cone(p);
1834: const index_type& pInd = this->_atlas->getIndex(patch, p);
1835: int j = -1;
1837: for(int i = pInd.prefix; i < pInd.prefix + pInd.index; ++i) {
1838: values[++j] = a[i];
1839: }
1840: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1841: const index_type& ind = this->_atlas->getIndex(patch, *p_iter);
1842: const int start = ind.prefix;
1843: const int length = ind.index;
1845: for(int i = start; i < start + length; ++i) {
1846: values[++j] = a[i];
1847: }
1848: }
1849: } else {
1850: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1851: int j = -1;
1853: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1854: const int start = i_iter->prefix;
1855: const int length = i_iter->index;
1857: for(int i = start; i < start + length; ++i) {
1858: values[++j] = a[i];
1859: }
1860: }
1861: }
1862: return values;
1863: };
1864: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {
1865: this->checkPatch(patch);
1866: // Using the index structure explicitly
1867: return &(this->_arrays[patch][this->_atlas->getIndex(patch, p).prefix]);
1868: };
1869: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
1870: this->checkPatch(patch);
1871: value_type *a = this->_arrays[patch];
1873: if (this->_atlas->getTopology()->height(patch) == 1) {
1874: // Only avoids the copy
1875: const Obj<typename sieve_type::coneSequence>& cone = this->_atlas->getTopology()->getPatch(patch)->cone(p);
1876: const index_type& pInd = this->_atlas->getIndex(patch, p);
1877: int j = -1;
1879: for(int i = pInd.prefix; i < pInd.prefix + pInd.index; ++i) {
1880: a[i] = v[++j];
1881: }
1882: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1883: const index_type& ind = this->_atlas->getIndex(patch, *p_iter);
1884: const int start = ind.prefix;
1885: const int length = ind.index;
1887: for(int i = start; i < start + length; ++i) {
1888: a[i] = v[++j];
1889: }
1890: }
1891: } else {
1892: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1893: int j = -1;
1895: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1896: const int start = i_iter->prefix;
1897: const int length = i_iter->index;
1899: for(int i = start; i < start + length; ++i) {
1900: a[i] = v[++j];
1901: }
1902: }
1903: }
1904: };
1905: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
1906: this->checkPatch(patch);
1907: value_type *a = this->_arrays[patch];
1909: if (this->_atlas->getTopology()->height(patch) == 1) {
1910: // Only avoids the copy
1911: const Obj<typename sieve_type::coneSequence>& cone = this->_atlas->getTopology()->getPatch(patch)->cone(p);
1912: const index_type& pInd = this->_atlas->getIndex(patch, p);
1913: int j = -1;
1915: for(int i = pInd.prefix; i < pInd.prefix + pInd.index; ++i) {
1916: a[i] += v[++j];
1917: }
1918: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1919: const index_type& ind = this->_atlas->getIndex(patch, *p_iter);
1920: const int start = ind.prefix;
1921: const int length = ind.index;
1923: for(int i = start; i < start + length; ++i) {
1924: a[i] += v[++j];
1925: }
1926: }
1927: } else {
1928: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1929: int j = -1;
1931: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1932: const int start = i_iter->prefix;
1933: const int length = i_iter->index;
1935: for(int i = start; i < start + length; ++i) {
1936: a[i] += v[++j];
1937: }
1938: }
1939: }
1940: };
1941: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
1942: this->checkPatch(patch);
1943: const index_type& ind = this->_atlas->getIndex(patch, p);
1944: value_type *a = &(this->_arrays[patch][ind.prefix]);
1946: // Using the index structure explicitly
1947: for(int i = 0; i < ind.index; ++i) {
1948: a[i] = v[i];
1949: }
1950: };
1951: template<typename Input>
1952: void update(const patch_type& patch, const point_type& p, const Obj<Input>& v) {
1953: this->checkPatch(patch);
1954: value_type *a = this->_arrays[patch];
1956: if (this->_atlas->getTopology()->height(patch) == 1) {
1957: // Only avoids the copy
1958: const Obj<typename sieve_type::coneSequence>& cone = this->_atlas->getTopology()->getPatch(patch)->cone(p);
1959: const index_type& pInd = this->_atlas->getIndex(patch, p);
1960: typename Input::iterator v_iter = v->begin();
1961: typename Input::iterator v_end = v->end();
1963: for(int i = pInd.prefix; i < pInd.prefix + pInd.index; ++i) {
1964: a[i] = *v_iter;
1965: ++v_iter;
1966: }
1967: for(typename sieve_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
1968: const index_type& ind = this->_atlas->getIndex(patch, *p_iter);
1969: const int start = ind.prefix;
1970: const int length = ind.index;
1972: for(int i = start; i < start + length; ++i) {
1973: a[i] = *v_iter;
1974: ++v_iter;
1975: }
1976: }
1977: } else {
1978: const Obj<typename atlas_type::IndexArray>& ind = this->_atlas->getIndices(patch, p);
1979: typename Input::iterator v_iter = v->begin();
1980: typename Input::iterator v_end = v->end();
1982: for(typename atlas_type::IndexArray::iterator i_iter = ind->begin(); i_iter != ind->end(); ++i_iter) {
1983: const int start = i_iter->prefix;
1984: const int length = i_iter->index;
1986: for(int i = start; i < start + length; ++i) {
1987: a[i] = *v_iter;
1988: ++v_iter;
1989: }
1990: }
1991: }
1992: };
1993: public: // Resizing
1994: void addPoint(const patch_type& patch, const point_type& point, const int size) {
1995: const typename atlas_type::chart_type& chart = this->_atlas->getChart(patch);
1997: if (chart.find(point) == chart.end()) {
1998: if (this->_atlasNew.isNull()) {
1999: this->_atlasNew = new atlas_type(this->_atlas->getTopology());
2000: this->_atlasNew->copy(this->_atlas);
2001: }
2002: this->_atlasNew->setFiberDimension(patch, point, size);
2003: }
2004: };
2005: void reallocate() {
2006: if (this->_atlasNew.isNull()) return;
2007: this->_atlasNew->orderPatches();
2008: const typename atlas_type::topology_type::sheaf_type& patches = this->_atlas->getTopology()->getPatches();
2010: for(typename atlas_type::topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2011: const typename atlas_type::patch_type& patch = p_iter->first;
2012: const typename atlas_type::chart_type& chart = this->_atlas->getChart(patch);
2013: const value_type *array = this->_arrays[patch];
2014: value_type *newArray = new value_type[this->_atlasNew->size(patch)];
2016: for(typename atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
2017: const int& size = c_iter->second.index;
2018: const int& offset = c_iter->second.prefix;
2019: const int& newOffset = this->_atlasNew->getIndex(patch, c_iter->first).prefix;
2021: for(int i = 0; i < size; ++i) {
2022: newArray[newOffset+i] = array[offset+i];
2023: }
2024: }
2025: delete [] this->_arrays[patch];
2026: this->_arrays[patch] = newArray;
2027: }
2028: this->_atlas = this->_atlasNew;
2029: this->_atlasNew = NULL;
2030: };
2031: public:
2032: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
2033: ostringstream txt;
2034: int rank;
2036: if (comm == MPI_COMM_NULL) {
2037: comm = this->comm();
2038: rank = this->commRank();
2039: } else {
2040: MPI_Comm_rank(comm, &rank);
2041: }
2042: if (name == "") {
2043: if(rank == 0) {
2044: txt << "viewing a Section" << std::endl;
2045: }
2046: } else {
2047: if(rank == 0) {
2048: txt << "viewing Section '" << name << "'" << std::endl;
2049: }
2050: }
2051: for(typename values_type::const_iterator a_iter = this->_arrays.begin(); a_iter != this->_arrays.end(); ++a_iter) {
2052: const patch_type patch = a_iter->first;
2053: const value_type *array = a_iter->second;
2055: txt << "[" << this->commRank() << "]: Patch " << patch << std::endl;
2056: const typename atlas_type::chart_type& chart = this->_atlas->getChart(patch);
2058: for(typename atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
2059: const typename atlas_type::index_type& idx = c_iter->second;
2061: if (idx.index != 0) {
2062: txt << "[" << this->commRank() << "]: " << c_iter->first << " dim " << idx.index << " offset " << idx.prefix << " ";
2063: for(int i = 0; i < idx.index; i++) {
2064: txt << " " << array[idx.prefix+i];
2065: }
2066: txt << std::endl;
2067: }
2068: }
2069: }
2070: PetscSynchronizedPrintf(comm, txt.str().c_str());
2071: PetscSynchronizedFlush(comm);
2072: };
2073: };
2075: template<typename Topology_, typename Value_>
2076: class ConstantSection : public ALE::ParallelObject {
2077: public:
2078: typedef Topology_ topology_type;
2079: typedef typename topology_type::patch_type patch_type;
2080: typedef typename topology_type::sieve_type sieve_type;
2081: typedef typename topology_type::point_type point_type;
2082: typedef Value_ value_type;
2083: protected:
2084: Obj<topology_type> _topology;
2085: const value_type _value;
2086: public:
2087: ConstantSection(MPI_Comm comm, const value_type value, const int debug = 0) : ParallelObject(comm, debug), _value(value) {
2088: this->_topology = new topology_type(comm, debug);
2089: };
2090: ConstantSection(const Obj<topology_type>& topology, const value_type value) : ParallelObject(topology->comm(), topology->debug()), _topology(topology), _value(value) {};
2091: virtual ~ConstantSection() {};
2092: public:
2093: void allocate() {};
2094: const value_type *restrict(const patch_type& patch) {return &this->_value;};
2095: // This should return something the size of the closure
2096: const value_type *restrict(const patch_type& patch, const point_type& p) {return &this->_value;};
2097: const value_type *restrictPoint(const patch_type& patch, const point_type& p) {return &this->_value;};
2098: void update(const patch_type& patch, const point_type& p, const value_type v[]) {
2099: throw ALE::Exception("Cannot update a ConstantSection");
2100: };
2101: void updateAdd(const patch_type& patch, const point_type& p, const value_type v[]) {
2102: throw ALE::Exception("Cannot update a ConstantSection");
2103: };
2104: void updatePoint(const patch_type& patch, const point_type& p, const value_type v[]) {
2105: throw ALE::Exception("Cannot update a ConstantSection");
2106: };
2107: template<typename Input>
2108: void update(const patch_type& patch, const point_type& p, const Obj<Input>& v) {
2109: throw ALE::Exception("Cannot update a ConstantSection");
2110: };
2111: public:
2112: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
2113: ostringstream txt;
2114: int rank;
2116: if (comm == MPI_COMM_NULL) {
2117: comm = this->comm();
2118: rank = this->commRank();
2119: } else {
2120: MPI_Comm_rank(comm, &rank);
2121: }
2122: if (name == "") {
2123: if(rank == 0) {
2124: txt << "viewing a ConstantSection with value " << this->_value << std::endl;
2125: }
2126: } else {
2127: if(rank == 0) {
2128: txt << "viewing ConstantSection '" << name << "' with value " << this->_value << std::endl;
2129: }
2130: }
2131: PetscSynchronizedPrintf(comm, txt.str().c_str());
2132: PetscSynchronizedFlush(comm);
2133: };
2134: };
2136: template<typename Point_>
2137: class DiscreteSieve {
2138: public:
2139: typedef Point_ point_type;
2140: typedef std::vector<point_type> coneSequence;
2141: typedef std::vector<point_type> coneSet;
2142: typedef std::vector<point_type> coneArray;
2143: typedef std::vector<point_type> supportSequence;
2144: typedef std::vector<point_type> supportSet;
2145: typedef std::vector<point_type> supportArray;
2146: typedef std::set<point_type> points_type;
2147: typedef points_type baseSequence;
2148: typedef points_type capSequence;
2149: protected:
2150: Obj<points_type> _points;
2151: Obj<coneSequence> _empty;
2152: Obj<coneSequence> _return;
2153: void _init() {
2154: this->_points = new points_type();
2155: this->_empty = new coneSequence();
2156: this->_return = new coneSequence();
2157: };
2158: public:
2159: DiscreteSieve() {
2160: this->_init();
2161: };
2162: template<typename Input>
2163: DiscreteSieve(const Obj<Input>& points) {
2164: this->_init();
2165: this->_points->insert(points->begin(), points->end());
2166: };
2167: virtual ~DiscreteSieve() {};
2168: public:
2169: void addPoint(const point_type& point) {
2170: this->_points->insert(point);
2171: };
2172: template<typename Input>
2173: void addPoints(const Obj<Input>& points) {
2174: this->_points->insert(points->begin(), points->end());
2175: };
2176: const Obj<coneSequence>& cone(const point_type& p) {return this->_empty;};
2177: template<typename Input>
2178: const Obj<coneSequence>& cone(const Input& p) {return this->_empty;};
2179: const Obj<coneSet>& nCone(const point_type& p, const int level) {
2180: if (level == 0) {
2181: return this->closure(p);
2182: } else {
2183: return this->_empty;
2184: }
2185: };
2186: const Obj<coneArray>& closure(const point_type& p) {
2187: this->_return->clear();
2188: this->_return->push_back(p);
2189: return this->_return;
2190: };
2191: const Obj<supportSequence>& support(const point_type& p) {return this->_empty;};
2192: template<typename Input>
2193: const Obj<supportSequence>& support(const Input& p) {return this->_empty;};
2194: const Obj<supportSet>& nSupport(const point_type& p, const int level) {
2195: if (level == 0) {
2196: return this->star(p);
2197: } else {
2198: return this->_empty;
2199: }
2200: };
2201: const Obj<supportArray>& star(const point_type& p) {
2202: this->_return->clear();
2203: this->_return->push_back(p);
2204: return this->_return;
2205: };
2206: const Obj<capSequence>& roots() {return this->_points;};
2207: const Obj<capSequence>& cap() {return this->_points;};
2208: const Obj<baseSequence>& leaves() {return this->_points;};
2209: const Obj<baseSequence>& base() {return this->_points;};
2210: template<typename Color>
2211: void addArrow(const point_type& p, const point_type& q, const Color& color) {
2212: throw ALE::Exception("Cannot add an arrow to a DiscreteSieve");
2213: };
2214: void stratify() {};
2215: void view(const std::string& name, MPI_Comm comm = MPI_COMM_NULL) const {
2216: ostringstream txt;
2217: int rank;
2219: if (comm == MPI_COMM_NULL) {
2220: comm = MPI_COMM_SELF;
2221: rank = 0;
2222: } else {
2223: MPI_Comm_rank(comm, &rank);
2224: }
2225: if (name == "") {
2226: if(rank == 0) {
2227: txt << "viewing a DiscreteSieve" << std::endl;
2228: }
2229: } else {
2230: if(rank == 0) {
2231: txt << "viewing DiscreteSieve '" << name << "'" << std::endl;
2232: }
2233: }
2234: for(typename points_type::const_iterator p_iter = this->_points->begin(); p_iter != this->_points->end(); ++p_iter) {
2235: txt << " Point " << *p_iter << std::endl;
2236: }
2237: PetscSynchronizedPrintf(comm, txt.str().c_str());
2238: PetscSynchronizedFlush(comm);
2239: };
2240: };
2243: template<typename Overlap_, typename Atlas_, typename Value_>
2244: class OverlapValues : public Section<Atlas_, Value_> {
2245: public:
2246: typedef typename Section<Atlas_, Value_>::patch_type patch_type;
2247: typedef typename Section<Atlas_, Value_>::topology_type topology_type;
2248: typedef Overlap_ overlap_type;
2249: typedef Atlas_ atlas_type;
2250: typedef Value_ value_type;
2251: typedef enum {SEND, RECEIVE} request_type;
2252: typedef std::map<patch_type, MPI_Request> requests_type;
2253: protected:
2254: int _tag;
2255: MPI_Datatype _datatype;
2256: requests_type _requests;
2257: public:
2258: OverlapValues(MPI_Comm comm, const int debug = 0) : Section<Atlas_, Value_>(comm, debug) {
2259: this->_tag = this->getNewTag();
2260: this->_datatype = this->getMPIDatatype();
2261: };
2262: OverlapValues(MPI_Comm comm, const int tag, const int debug) : Section<Atlas_, Value_>(comm, debug), _tag(tag) {
2263: this->_datatype = this->getMPIDatatype();
2264: };
2265: virtual ~OverlapValues() {};
2266: protected:
2267: MPI_Datatype getMPIDatatype() {
2268: if (sizeof(value_type) == 4) {
2269: return MPI_INT;
2270: } else if (sizeof(value_type) == 8) {
2271: return MPI_DOUBLE;
2272: } else if (sizeof(value_type) == 28) {
2273: int blen[2];
2274: MPI_Aint indices[2];
2275: MPI_Datatype oldtypes[2], newtype;
2276: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_INT;
2277: blen[1] = 3; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
2278: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
2279: MPI_Type_commit(&newtype);
2280: return newtype;
2281: }
2282: throw ALE::Exception("Cannot determine MPI type for value type");
2283: };
2284: int getNewTag() {
2285: static int tagKeyval = MPI_KEYVAL_INVALID;
2286: int *tagvalp = NULL, *maxval, flg;
2288: if (tagKeyval == MPI_KEYVAL_INVALID) {
2289: PetscMalloc(sizeof(int), &tagvalp);
2290: MPI_Keyval_create(MPI_NULL_COPY_FN, Petsc_DelTag, &tagKeyval, (void *) NULL);
2291: MPI_Attr_put(this->_comm, tagKeyval, tagvalp);
2292: tagvalp[0] = 0;
2293: }
2294: MPI_Attr_get(this->_comm, tagKeyval, (void **) &tagvalp, &flg);
2295: if (tagvalp[0] < 1) {
2296: MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, (void **) &maxval, &flg);
2297: tagvalp[0] = *maxval - 128; // hope that any still active tags were issued right at the beginning of the run
2298: }
2299: //std::cout << "[" << this->commRank() << "]Got new tag " << tagvalp[0] << std::endl;
2300: return tagvalp[0]--;
2301: };
2302: public: // Accessors
2303: int getTag() const {return this->_tag;};
2304: void setTag(const int tag) {this->_tag = tag;};
2305: public:
2306: void construct(const int size) {
2307: const typename topology_type::sheaf_type& patches = this->_atlas->getTopology()->getPatches();
2309: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2310: const Obj<typename topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
2311: int rank = p_iter->first;
2313: for(typename topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
2314: this->_atlas->setFiberDimension(rank, *b_iter, size);
2315: }
2316: }
2317: };
2318: template<typename Sizer>
2319: void construct(const Obj<Sizer>& sizer) {
2320: const typename topology_type::sheaf_type& patches = this->_atlas->getTopology()->getPatches();
2322: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2323: const Obj<typename topology_type::sieve_type::baseSequence>& base = p_iter->second->base();
2324: int rank = p_iter->first;
2326: for(typename topology_type::sieve_type::baseSequence::iterator b_iter = base->begin(); b_iter != base->end(); ++b_iter) {
2327: this->_atlas->setFiberDimension(rank, *b_iter, *(sizer->restrict(rank, *b_iter)));
2328: }
2329: }
2330: };
2331: void constructCommunication(const request_type& requestType) {
2332: const typename topology_type::sheaf_type& patches = this->getAtlas()->getTopology()->getPatches();
2334: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2335: const patch_type patch = p_iter->first;
2336: MPI_Request request;
2338: if (requestType == RECEIVE) {
2339: if (this->_debug) {std::cout <<"["<<this->commRank()<<"] Receiving data(" << this->_atlas->size(patch) << ") from " << patch << " tag " << this->_tag << std::endl;}
2340: MPI_Recv_init(this->_arrays[patch], this->_atlas->size(patch), this->_datatype, patch, this->_tag, this->_comm, &request);
2341: } else {
2342: if (this->_debug) {std::cout <<"["<<this->commRank()<<"] Sending data (" << this->_atlas->size(patch) << ") to " << patch << " tag " << this->_tag << std::endl;}
2343: MPI_Send_init(this->_arrays[patch], this->_atlas->size(patch), this->_datatype, patch, this->_tag, this->_comm, &request);
2344: }
2345: this->_requests[patch] = request;
2346: }
2347: };
2348: void startCommunication() {
2349: const typename topology_type::sheaf_type& patches = this->getAtlas()->getTopology()->getPatches();
2351: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2352: MPI_Request request = this->_requests[p_iter->first];
2354: MPI_Start(&request);
2355: }
2356: };
2357: void endCommunication() {
2358: const typename topology_type::sheaf_type& patches = this->getAtlas()->getTopology()->getPatches();
2359: MPI_Status status;
2361: for(typename topology_type::sheaf_type::const_iterator p_iter = patches.begin(); p_iter != patches.end(); ++p_iter) {
2362: MPI_Request request = this->_requests[p_iter->first];
2364: MPI_Wait(&request, &status);
2365: }
2366: };
2367: };
2368: }
2369: }
2371: #endif