Actual source code: Hierarchy_New.hh
1: /*
2: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
4: Hierarchy_New.hh:
6: Routines for coarsening arbitrary simplicial meshes expressed as sieves using
7: a localized modification of the algorithm in Miller (1999)
9: - Peter Brune
11: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
12: */
14: #include <list>
15: #include <stdlib.h>
17: #include <Mesh.hh>
18: #include "SieveAlgorithms.hh"
19: #include "MeshSurgery.hh"
20: #include <SieveAlgorithms.hh>
21: #include <Selection.hh>
27: /*
28: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
29: Hierarchy_createCoarsenVertexSet:
31: Inputs:
33: original_mesh: the mesh in question
34: spacing: the spacing function over the vertices of the original_mesh
35: interior_set: the points that could potentially be pruned
36: boundary_set: the points that are forced into the mesh
37: beta: the expansion of the spacing function
38: maxIndex: the maximum index in original_mesh, used for avoiding tag collisions
40: Returns:
42: The the set of vertices that are left standing after the coarsening procedure
44: Remarks:
46: Doesn't modify the original_mesh topology, works only with interactions between
47: the boundary set and the interior set and within the interior set. Use for all
48: coarsening and then feed the output into a remesher. Can be used for each
49: layer of boundary coarsening for 2D/3D/WhateverD meshes.
51: Builds a rather expensive connection_sieve structure in the process. This
52: has no labels associated with it, so it's less expensive than the mesh is.
53: Subsets of the mesh could be coarsened at a time with an internal skeleton
54: boundary if this is too expensive.
57: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
58: */
60: ALE::Obj<ALE::Mesh::sieve_type::supportSet> Hierarchy_CoarsenVertexSet(ALE::Obj<ALE::Mesh> original_mesh,
61: ALE::Obj<ALE::Mesh::real_section_type> spacing,
62: ALE::Obj<ALE::Mesh::sieve_type::supportSet> interior_set,
63: ALE::Obj<ALE::Mesh::sieve_type::supportSet> boundary_set,
64: double beta,
65: ALE::Mesh::point_type maxIndex,
66: int * comparisons = PETSC_NULL) {
67:
70: //PetscPrintf(original_mesh->comm(), "%f\n", beta);
71: ALE::Obj<ALE::Mesh::sieve_type> original_sieve = original_mesh->getSieve();
72: ALE::Obj<ALE::Mesh::real_section_type> coordinates = original_mesh->getRealSection("coordinates");
74: ALE::Mesh::point_type cur_new_index = maxIndex+1;
76: //build the coarsening structure with the following links: links within the interior_set and links between the boundary_set and interior_set
78: ALE::Obj<ALE::Mesh::sieve_type> connection_sieve = new ALE::Mesh::sieve_type(original_sieve->comm(), original_sieve->debug());
79: ALE::Obj<ALE::Mesh> connection_mesh = new ALE::Mesh(original_mesh->comm(), original_mesh->debug());
80: connection_mesh->setSieve(connection_sieve);
83: ALE::Mesh::sieve_type::supportSet::iterator is_iter = interior_set->begin();
84: ALE::Mesh::sieve_type::supportSet::iterator is_iter_end = interior_set->end();
85: ALE::Mesh::sieve_type::supportSet::iterator bs_iter_end = boundary_set->end();
87: ALE::Obj<ALE::Mesh::sieve_type::supportSet> output_vertices = new ALE::Mesh::sieve_type::supportSet();
89: std::list<ALE::Mesh::point_type> comparison_queue; //edges go here
91: int dim = original_mesh->getDimension();
93: double *a_coords = new double[dim], *b_coords = new double[dim];
95: ALE::Mesh::sieve_type::supportSet candidate_vertices;
97: ALE::Obj<ALE::Mesh::sieve_type::supportSet> link;
98: ALE::Obj<ALE::Mesh::sieve_type::supportSet> line = ALE::Mesh::sieve_type::supportSet();
100: while (is_iter != is_iter_end) {
101: //create the set of candidate points
102: if (boundary_set->find(*is_iter) == bs_iter_end) {
103: candidate_vertices.insert(*is_iter);
104: } else {
105: output_vertices->insert(*is_iter);
106: }
108: ALE::Obj<ALE::Mesh::sieve_type::coneSet> neighbors = original_sieve->cone(original_sieve->support(*is_iter));
109: ALE::Mesh::sieve_type::coneSet::iterator n_iter = neighbors->begin();
110: ALE::Mesh::sieve_type::coneSet::iterator n_iter_end = neighbors->end();
111: while (n_iter != n_iter_end) {
112: if (*n_iter != *is_iter)
113: if (boundary_set->find(*n_iter) != bs_iter_end) {
114: connection_sieve->addArrow(*n_iter, cur_new_index);
115: connection_sieve->addArrow(*is_iter, cur_new_index);
116: comparison_queue.push_back(cur_new_index);
117: cur_new_index++;
118: } else if (interior_set->find(*n_iter) != interior_set->end()) {
119: if (candidate_vertices.find(*n_iter) == candidate_vertices.end()) {
120: connection_sieve->addArrow(*n_iter, cur_new_index);
121: connection_sieve->addArrow(*is_iter, cur_new_index);
122: cur_new_index++;
123: }
124: }
125: n_iter++;
126: }
127: is_iter++;
128: }
129: // PetscPrintf(original_mesh->comm(), "current new index: %d\n", cur_new_index);
130: //PetscPrintf(original_mesh->comm(), "trying to stratify\n");
131: connection_mesh->stratify();
132: output_vertices->clear();
133: //PetscPrintf(original_mesh->comm(), "onto decimation!\n");
134: //decimate the new structure
135: while (!candidate_vertices.empty()) {
136: //PetscPrintf(original_mesh->comm(), "starting with this round.\n");
137: while (!comparison_queue.empty()) {
138: ALE::Mesh::point_type current_edge = comparison_queue.front();
139: comparison_queue.pop_front();
140: //logic for comparisons
141: ALE::Obj<ALE::Mesh::sieve_type::coneSequence> endpoints = connection_sieve->cone(current_edge);
142: //PetscPrintf(original_mesh->comm(), "got the endpoints: size %d\n", endpoints->size());
143: // if (endpoints->size() != 2) throw ALE::Exception("not the connection sieve?");
144: ALE::Mesh::point_type a, b;
145: bool a_included = false, b_included = false;
146: if (endpoints->size() == 2) {
147: ALE::Mesh::sieve_type::coneSequence::iterator ep_iter = endpoints->begin();
148: a = *ep_iter;
149: ep_iter++;
150: b = *ep_iter;
151: if (output_vertices->find(a) != output_vertices->end()) {
152: a_included = true;
153: //PetscPrintf(original_mesh->comm(), "%d is in output\n", a);
154: } else if (boundary_set->find(a) != boundary_set->end()) {
155: a_included = true;
156: //PetscPrintf(original_mesh->comm(), "%d is in boundary\n ", a);
157: }
158: if (output_vertices->find(b) != output_vertices->end()) {
159: //PetscPrintf(original_mesh->comm(), "%d is in output\n",b);
160: b_included = true;
161: } else if (boundary_set->find(b) != boundary_set->end()) {
162: //PetscPrintf(original_mesh->comm(), "%d is in boundary\n",b);
163: b_included = true;
164: }
165: }
166: if ((a_included && b_included) || ((!b_included) && (!a_included))) {
167: //PetscPrintf(original_mesh->comm(), "edge is irrelevant\n");
168: //either both are included or both are not included so this comparison doesn't need to be done
169: } else {
170: //a should be the point that is included in the mesh already; if it is not; swap it
171: if (b_included) {
172: ALE::Mesh::point_type swap_vertex = a;
173: a = b;
174: b = swap_vertex;
175: }
176: //do the comparison
177: PetscMemcpy(a_coords, coordinates->restrictPoint(a), dim*sizeof(double));
178: PetscMemcpy(b_coords, coordinates->restrictPoint(b), dim*sizeof(double));
179: double space_a = *spacing->restrictPoint(a);
180: double space_b = *spacing->restrictPoint(b);
181: double space = beta*(space_a + space_b);
182: double dist = 0;
183: for (int i = 0; i < dim; i++) {
184: dist += (a_coords[i] - b_coords[i])*(a_coords[i] - b_coords[i]);
185: }
186: dist = sqrt(dist);
187: if (dist < space) { //collision case; remove b from the candidate points and the sieve.
188: //PetscPrintf(original_mesh->comm(), "removal criterion satisfied\n");
189: ALE::Mesh::sieve_type::supportSet::iterator removal_iter = candidate_vertices.find(b);
190: if (removal_iter != candidate_vertices.end()) {
191: //remove b from the queue, reconnect link vertices to a
192: ALE::Obj<ALE::Mesh::sieve_type::supportSequence> b_edges = connection_sieve->support(b);
193: ALE::Mesh::sieve_type::supportSequence::iterator be_iter = b_edges->begin();
194: ALE::Mesh::sieve_type::supportSequence::iterator be_iter_end = b_edges->end();
195: ALE::Mesh::sieve_type::supportSet remove_these_points;
196: while (be_iter != be_iter_end) {
197: connection_sieve->addArrow(a, *be_iter);
198: if (connection_sieve->cone(*be_iter)->size() != 3) { //SHOULD contain a, b, and the associated link vertex; can be a and b, or just
199: remove_these_points.insert(*be_iter);
200: //connection_sieve->removeBasePoint(*be_iter);
201: //connection_sieve->removeCapPoint(*be_iter);
202:
203: } else { //readd the edges to the queue
204: comparison_queue.push_back(*be_iter);
205: }
206: be_iter++;
207: }
208: for (ALE::Mesh::sieve_type::supportSet::iterator rtp_iter = remove_these_points.begin(); rtp_iter != remove_these_points.end(); rtp_iter++) {
209: connection_sieve->removeBasePoint(*rtp_iter);
210: connection_sieve->removeCapPoint(*rtp_iter);
211: }
212: connection_sieve->removeBasePoint(b);
213: connection_sieve->removeCapPoint(b);
214: //remove it from the candidate points
215: candidate_vertices.erase(removal_iter);
216: //PetscPrintf(original_mesh->comm(), "removed a vertex\n");
217: } //end of removal_iter if
218: } //end of spacing if
219: } //end of else denoting that this comparison matters
220: } //end of !comparison_queue.empty() while
221: //the comparison queue is empty; the present output_vertices is compatible
222: //take one point from the candidate_vertices and move it to the output_vertices, adding its edges to the comparison queue
223: if (!candidate_vertices.empty()) { //we could have removed the final vertex in the previous coarsening. check
224: ALE::Mesh::sieve_type::supportSet::iterator new_added_vertex_iter = candidate_vertices.begin();
225: ALE::Mesh::point_type new_added_vertex = *new_added_vertex_iter;
226: candidate_vertices.erase(new_added_vertex_iter);
227: output_vertices->insert(new_added_vertex); //add this new vertex to the set of vertices in the output.
228: ALE::Obj<ALE::Mesh::sieve_type::supportSequence> nav_support = connection_sieve->support(new_added_vertex);
229: ALE::Mesh::sieve_type::supportSequence::iterator nav_iter = nav_support->begin();
230: ALE::Mesh::sieve_type::supportSequence::iterator nav_iter_end = nav_support->end();
231: while (nav_iter != nav_iter_end) {
232: comparison_queue.push_back(*nav_iter);
233: nav_iter++;
234: }
235: } //end of candidate_vertices.empty() if for adding a new vertex
236: } //end of !candidate_vertices.empty() while
237: //PetscPrintf(original_mesh->comm(), "%d vertices included in coarse set\n", output_vertices->size());
238: delete [] a_coords; delete [] b_coords;
239: return output_vertices;
240: }
241: /*
242: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
243: Hierarchy_createEffective2DBoundary:
245: Use:
247: Returns a noninterpolated 2D boundary mesh, including faults, for use in coarsening
249: Inputs:
251: original_mesh: an effectively 3D simplicial mesh
252: maxIndex: the beginning of the free indices for use in the boundary mesh
253: forced_bound_mesh: an optional internal fault mesh (uninterpolated)
255: Output:
257: An uninterpolated boundary mesh containing all exterior faces in the 3D mesh, as well
258: as the additional internal fault mesh given by forced_bound_mesh
260: Assumptions:
262: forced_bound_mesh is uninterpolated and distinct from the exterior
264: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
265: */
269: ALE::Obj<ALE::Mesh> Hierarchy_createEffective2DBoundary (ALE::Obj<ALE::Mesh> original_mesh, ALE::Mesh::point_type maxIndex, ALE::Obj<ALE::Mesh> forced_bound_mesh = PETSC_NULL) {
270: //create a set of "essential" faces
271: //from a given 3D mesh
272: int depth = original_mesh->depth();
273: int dim = original_mesh->getDimension();
274: int cur_available_index = maxIndex+1;
275: ALE::Obj<ALE::Mesh::sieve_type::supportSet> line = new ALE::Mesh::sieve_type::supportSet();
276: ALE::Obj<ALE::Mesh::sieve_type::supportSet> link;
277: ALE::Obj<ALE::Mesh::sieve_type> original_sieve = original_mesh->getSieve();
279: //setup the output
281: //build the exterior boundary
282: ALE::Obj<ALE::Mesh> output_mesh = new ALE::Mesh(original_mesh->comm(), dim, original_mesh->debug());
283: ALE::Obj<ALE::Mesh::sieve_type> output_sieve = new ALE::Mesh::sieve_type(original_mesh->comm(), original_mesh->debug());
284: output_mesh->setSieve(output_sieve);
286: if (depth == 1) { //noninterpolated case; for each cell look at the associated vertex subsets and assume simplicial
287: //PetscPrintf(original_mesh->comm(), "Uninterpolated case\n");
288: ALE::Obj<ALE::Mesh::label_sequence> cells = original_mesh->heightStratum(0);
289: ALE::Mesh::label_sequence::iterator c_iter = cells->begin();
290: ALE::Mesh::label_sequence::iterator c_iter_end = cells->end();
291: while (c_iter != c_iter_end) {
292: ALE::Obj<ALE::Mesh::sieve_type::coneSequence> cell_corners = original_sieve->cone(*c_iter);
293: ALE::Mesh::sieve_type::coneSequence::iterator cc_iter = cell_corners->begin();
294: ALE::Mesh::sieve_type::coneSequence::iterator cc_iter_end = cell_corners->end();
295: line->clear();
296: while (cc_iter != cc_iter_end) {
297: line->insert(*cc_iter);
298: cc_iter++;
299: }
300: //PetscPrintf(original_mesh->comm(), "%d vertices in the cone\n", cell_corners->size());
301: //assumed simplicial; each n-1 subset is a valid face
302: cc_iter = cell_corners->begin();
303: while (cc_iter != cc_iter_end) {
304: line->erase(line->find(*cc_iter));
305: link = original_sieve->nJoin1(line);
306: // PetscPrintf(original_mesh->comm(), "%d vertices in the face\n", line->size());
307: if (line->size() != 3) throw ALE::Exception("bad line");
308: if (link->size() == 1) {
309: for (ALE::Mesh::sieve_type::supportSet::iterator f_iter = line->begin(); f_iter != line->end(); f_iter++) {
310: output_sieve->addArrow(*f_iter, cur_available_index);
311: }
312: cur_available_index++;
313: }
314: line->insert(*cc_iter);
315: cc_iter++;
316: }
317: c_iter++;
318: }
319: } else { //interpolated case; we have the faces, just add them
320: //PetscPrintf(original_mesh->comm(), "Interpolated case\n");
321: ALE::Obj<ALE::Mesh::label_sequence> outward_faces = original_mesh->heightStratum(1); //assume wlog the faces are heightstratum 1
322: ALE::Mesh::label_sequence::iterator f_iter = outward_faces->begin();
323: ALE::Mesh::label_sequence::iterator f_iter_end = outward_faces->end();
324: while(f_iter != f_iter_end) {
325: //PetscPrintf(original_mesh->comm(), "testing...\n");
326: if (original_sieve->support(*f_iter)->size() != 2) {//add the uninterpolated face to the boundary.
327: //if (cur_available_index <= *f_iter) cur_available-index = *f_iter+1;
328: ALE::Obj<ALE::Mesh::sieve_type::coneArray> corners = original_sieve->nCone(*f_iter, depth-1);
329: ALE::Mesh::sieve_type::coneArray::iterator c_iter = corners->begin();
330: ALE::Mesh::sieve_type::coneArray::iterator c_iter_end = corners->end();
331: while (c_iter != c_iter_end) {
332: output_sieve->addArrow(*c_iter, cur_available_index);
333: c_iter++;
334: }
335: cur_available_index++;
336: }
337: f_iter++;
338: }
339: }
340: //PetscPrintf(original_mesh->comm(), "Done with finding the boundary\n");
341: //force in the forced_bound_mesh, which will include any included faults
343: if (forced_bound_mesh) {
344: if (forced_bound_mesh->depth() != 1) throw ALE::Exception("Needs noninterpolated boundary meshes");
345: ALE::Obj<ALE::Mesh::label_sequence> forced_boundary_cells = forced_bound_mesh->heightStratum(0);
346: ALE::Mesh::label_sequence::iterator fbc_iter = forced_boundary_cells->begin();
347: ALE::Mesh::label_sequence::iterator fbc_iter_end = forced_boundary_cells->end();
348: while (fbc_iter != fbc_iter_end) {
349: // if (*fbc_iter > cur_available_index) cur_available_index = *fbc_iter+1;
350: ALE::Obj<ALE::Mesh::sieve_type::coneSequence> corners = forced_bound_mesh->getSieve()->cone(*fbc_iter);
351: ALE::Mesh::sieve_type::coneSequence::iterator c_iter = corners->begin();
352: ALE::Mesh::sieve_type::coneSequence::iterator c_iter_end = corners->end();
353: while (c_iter != c_iter_end) {
354: output_sieve->addArrow(*c_iter, cur_available_index);
355: c_iter++;
356: }
357: cur_available_index++;
358: fbc_iter++;
359: }
360: }
361: output_mesh->stratify();
362: output_mesh->setRealSection("coordinates", original_mesh->getRealSection("coordinates"));
363: //PetscPrintf(output_mesh->comm(), "done with 2D boundary find.\n");
364: return output_mesh;
365: }
367: /*
368: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
369: Hierarchy_createEffective1DBoundary:
371: Use:
373: Returns a noninterpolated 1D boundary mesh, including faults, for use in coarsening
375: Inputs:
377: original_mesh: an effectively 2D simplicial mesh
378: maxIndex: the beginning of the free indices for use in the boundary mesh
379: forced_bound_mesh: an optional internal fault mesh (uninterpolated)
381: Output:
383: An uninterpolated boundary mesh containing all exterior edges in the 2D mesh, as well
384: as the additional internal fault mesh given by forced_bound_mesh
386: Assumptions:
388: forced_bound_mesh is uninterpolated and distinct from the exterior
390: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
391: */
395: ALE::Obj<ALE::Mesh> Hierarchy_createEffective1DBoundary (ALE::Obj<ALE::Mesh> original_mesh, ALE::Mesh::point_type maxIndex, ALE::Obj<ALE::Mesh> forced_bound_mesh = PETSC_NULL) {
396: //create a set of "essential" edges
397: //from a given 2D mesh or 3D boundary mesh
398: int dim = original_mesh->getDimension();
399: ALE::Obj<ALE::Mesh::sieve_type> original_sieve = original_mesh->getSieve();
400: ALE::Obj<ALE::Mesh::real_section_type> coordinates = original_mesh->getRealSection("coordinates");
401: double *a_coords = new double[dim], *b_coords = new double[dim], *c_coords = new double[dim], *d_coords = new double[dim];
402: int depth = original_mesh->depth();
403: ALE::Mesh::point_type cur_available_index = maxIndex+1;
404: ALE::Obj<ALE::Mesh> output_mesh = new ALE::Mesh(original_mesh->comm(), dim, original_mesh->debug());
405: ALE::Obj<ALE::Mesh::sieve_type> output_sieve = new ALE::Mesh::sieve_type(original_mesh->comm(), original_mesh->debug());
406: output_mesh->setSieve(output_sieve);
407: ALE::Obj<ALE::Mesh::label_type> marker = output_mesh->createLabel("marker"); //EVERYTHING here gets the marker name
408: ALE::Obj<ALE::Mesh::sieve_type::supportSet> line = new ALE::Mesh::sieve_type::supportSet();
409: ALE::Obj<ALE::Mesh::sieve_type::supportSet> link;
411: //force in the forced boundary mesh (2D)
412: if (forced_bound_mesh) {
413: //PetscPrintf(original_mesh->comm(), "forcing in the boundary mesh\n");
414: if (forced_bound_mesh->depth() != 1) throw ALE::Exception("Needs noninterpolated boundary mesh"); //haha doesn't matter, but still
415: ALE::Obj<ALE::Mesh::label_sequence> bound_edges = forced_bound_mesh->heightStratum(0);
416: ALE::Mesh::label_sequence::iterator be_iter = bound_edges->begin();
417: ALE::Mesh::label_sequence::iterator be_iter_end = bound_edges->end();
418: while (be_iter != be_iter_end) {
419: ALE::Obj<ALE::Mesh::sieve_type::coneSequence> edge_ends = forced_bound_mesh->getSieve()->cone(*be_iter);
420: ALE::Mesh::sieve_type::coneSequence::iterator ee_iter = edge_ends->begin();
421: ALE::Mesh::sieve_type::coneSequence::iterator ee_iter_end = edge_ends->end();
422: while (ee_iter != ee_iter_end) {
423: output_sieve->addArrow(*ee_iter, cur_available_index);
424: ee_iter++;
425: }
426: cur_available_index++;
427: be_iter++;
428: }
429: }
430: //PetscPrintf(original_mesh->comm(), "Forcing in edges\n");
431: //find topologically needed edges -- interpolated means support != 2; noninterpolated means njoin1 != 2
432: //find geometrically attractive edges (3D); doublet-normals above a certain threshhold.
433: if (depth == 1) { //uninterpolated case
434: //PetscPrintf(original_mesh->comm(), "Uninterpolated\n");
435: ALE::Mesh::sieve_type::supportSet seen;
436: seen.clear();
437: ALE::Obj<ALE::Mesh::label_sequence> vertices = original_mesh->depthStratum(0);
438: ALE::Mesh::label_sequence::iterator v_iter = vertices->begin();
439: ALE::Mesh::label_sequence::iterator v_iter_end = vertices->end();
440: while (v_iter != v_iter_end) {
441: seen.insert(*v_iter);
442: ALE::Obj<ALE::Mesh::sieve_type::coneSet> neighbors = original_sieve->cone(original_sieve->support(*v_iter));
443: ALE::Mesh::sieve_type::coneSet::iterator n_iter = neighbors->begin();
444: ALE::Mesh::sieve_type::coneSet::iterator n_iter_end = neighbors->end();
445: while (n_iter != n_iter_end) { if (*n_iter != *v_iter) {
446: if (seen.find(*n_iter) == seen.end()) { //if we've seen this vertex and therefore done this edge before..
447: line->clear();
448: line->insert(*n_iter);
449: line->insert(*v_iter);
450: link = original_sieve->nJoin1(line);
451: if (link->size() != 2) {
452: //either a boundary in 2D or a fault/mesh intersection in 3D
453: output_sieve->addArrow(*v_iter, cur_available_index);
454: output_sieve->addArrow(*n_iter, cur_available_index);
455: cur_available_index++;
456: }else if (dim == 3) { //do the extra tests for doublet-wide principal curvature
457: //build the doublet
458: ALE::Mesh::point_type a, b, c, d;
459: a = *n_iter;
460: b = *v_iter;
461: ALE::Obj<ALE::Mesh::sieve_type::coneSet> doublet_corners = original_sieve->cone(link);
462: ALE::Mesh::sieve_type::coneSet::iterator dc_iter = doublet_corners->begin();
463: ALE::Mesh::sieve_type::coneSet::iterator dc_iter_end = doublet_corners->end();
464: int corner_number = 0;
465: while (dc_iter != dc_iter_end) {
466: if (*dc_iter != a && *dc_iter != b) {
467: if (corner_number == 0) {
468: c = *dc_iter;
469: corner_number = 1;
470: } else {
471: d = *dc_iter;
472: }
473: }
474: dc_iter++;
475: }
476: PetscMemcpy(a_coords, coordinates->restrictPoint(a), dim*sizeof(double));
477: PetscMemcpy(b_coords, coordinates->restrictPoint(b), dim*sizeof(double));
478: PetscMemcpy(c_coords, coordinates->restrictPoint(c), dim*sizeof(double));
479: PetscMemcpy(d_coords, coordinates->restrictPoint(d), dim*sizeof(double));
480: double d_angle = ALE::doublet_angle(dim, a_coords, b_coords, c_coords, d_coords);
481: //PetscPrintf(original_mesh->comm(), "doublet angle: %f\n", d_angle);
482: if (d_angle > 0.79) { //arbitrary cutoff of around 45 degrees
483: output_sieve->addArrow(*v_iter, cur_available_index);
484: output_sieve->addArrow(*n_iter, cur_available_index);
485: cur_available_index++;
486: }
487: }
488: }
489: }
490: n_iter++;
491: }
492: v_iter++;
493: }
494: } else { //interpolated case
495: //PetscPrintf(original_mesh->comm(), "Interpolated\n");
496: ALE::Obj<ALE::Mesh::label_sequence> edges = original_mesh->depthStratum(1);
497: ALE::Mesh::label_sequence::iterator e_iter = edges->begin();
498: ALE::Mesh::label_sequence::iterator e_iter_end = edges->end();
499: while (e_iter != e_iter_end) {
500: //PetscPrintf(original_mesh->comm(), "in the loop...\n");
501: ALE::Obj<ALE::Mesh::sieve_type::supportSequence> edge_support = original_sieve->support(*e_iter);
502: if (edge_support->size() != 2) {
503: ALE::Obj<ALE::Mesh::sieve_type::coneSequence> edge_cone = original_sieve->cone(*e_iter);
504: ALE::Mesh::sieve_type::coneSequence::iterator ec_iter = edge_cone->begin();
505: ALE::Mesh::sieve_type::coneSequence::iterator ec_iter_end = edge_cone->end();
506: while (ec_iter != ec_iter_end) {
507: output_sieve->addArrow(*ec_iter, cur_available_index);
508: ec_iter++;
509: }
510: cur_available_index++;
511: } else if (dim == 3) { //boundary meshes will be UNINTERPOLATED, so this won't come up
512: }
513: e_iter++;
514: }
515: }
516: output_mesh->stratify();
517: ALE::Obj<ALE::Mesh::label_sequence> points = output_mesh->depthStratum(0); //boundary vertices;
518: ALE::Mesh::label_sequence::iterator p_iter = points->begin();
519: ALE::Mesh::label_sequence::iterator p_iter_end = points->end();
520: while (p_iter != p_iter_end) {
521: output_mesh->setValue(marker, *p_iter, 1);
522: p_iter++;
523: }
524: points = output_mesh->heightStratum(0); //boundary edges;
525: p_iter = points->begin();
526: p_iter_end = points->end();
527: while (p_iter != p_iter_end) {
528: output_mesh->setValue(marker, *p_iter, 1);
529: p_iter++;
530: }
531: output_mesh->setRealSection("coordinates", original_mesh->getRealSection("coordinates"));
532: output_mesh->stratify();
533: //PetscPrintf(output_mesh->comm(), "leaving 1D Boundary Building\n");
534: delete [] a_coords; delete [] b_coords; delete [] c_coords; delete [] d_coords;
535: return output_mesh;
536: }
538: /*
539: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
540: Hierarchy_createEffective0DBoundary:
542: Use:
544: Returns a set of topologically necessary vertices from an effectively 1D mesh
546: Inputs:
548: original_mesh: an effectively 1D simplicial mesh
550: Output:
552: A set of vertices from the 1D simplicial mesh that are topologically necessary
555: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
556: */
561: ALE::Obj<ALE::Mesh::sieve_type::supportSet> Hierarchy_createEffective0DBoundary(ALE::Obj<ALE::Mesh> original_mesh) {
562: ALE::Obj<ALE::Mesh::real_section_type> coordinates = original_mesh->getRealSection("coordinates");
563: //PetscPrintf(original_mesh->comm(), "In 0D Boundary Building\n");
564: int dim = original_mesh->getDimension();
565: double *a_coords = new double[dim], *b_coords = new double[dim], *c_coords = new double[dim];
566: //create a set of "essential" vertices from the 1D boundary meshes given
567: //find the topologically necessary vertices in the 1D mesh -- this means vertices on which an edge terminates or is noncontractible.
568: ALE::Obj<ALE::Mesh::sieve_type::supportSet> output_set = new ALE::Mesh::sieve_type::supportSet();
569: ALE::Obj<ALE::Mesh::sieve_type> original_sieve = original_mesh->getSieve();
570: //go through the vertices of this set looking for ones with support size greater than 3.
571: ALE::Obj<ALE::Mesh::label_sequence> vertices = original_mesh->depthStratum(0);
572: ALE::Mesh::label_sequence::iterator v_iter = vertices->begin();
573: ALE::Mesh::label_sequence::iterator v_iter_end = vertices->end();
574: while (v_iter != v_iter_end) {
575: if (original_sieve->support(*v_iter)->size() != 2) {
576: output_set->insert(*v_iter);
577: } else if (dim > 1) { //angles
578: //if the angle between eges is greater than say... 45 degrees, include it.
579: ALE::Obj<ALE::Mesh::sieve_type::coneSet> neighbors = original_sieve->cone(original_sieve->support(*v_iter));
580: ALE::Mesh::sieve_type::coneSet::iterator n_iter = neighbors->begin();
581: ALE::Mesh::point_type b, c;
582: //should have two neighbors
583: if (*n_iter != *v_iter) b = *n_iter;
584: n_iter++;
585: if (*n_iter != *v_iter) { c = b; b = *n_iter;}
586: n_iter++;
587: if (*n_iter != *v_iter) { c = b; b = *n_iter;}
588: PetscMemcpy(a_coords, coordinates->restrictPoint(*v_iter), dim*sizeof(double));
589: PetscMemcpy(b_coords, coordinates->restrictPoint(b), dim*sizeof(double));
590: PetscMemcpy(c_coords, coordinates->restrictPoint(c), dim*sizeof(double));
591: double angle = ALE::corner_angle(dim, a_coords, b_coords, c_coords);
592: if (fabs(angle - 3.14159) > 0.78) output_set->insert(*v_iter);
593: }
594: v_iter++;
595: }
596: //we can't do curvatures here; leave this to another thing we'll include in this set later
597: //PetscPrintf(original_mesh->comm(), "leaving 0D Boundary Building\n");
598: delete [] a_coords; delete [] b_coords; delete [] c_coords;
599: return output_set;
600: }
602: /*
603: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
604: Hierarchy_curvatureSet:
606: Use:
608: Returns a set of geometrically attractive vertices from an effectively 1D or 2D mesh
610: Inputs:
612: original_mesh: an effectively 2D simplicial mesh
615: Output:
617: A set of vertices from the 2D simplicial mesh that are geometrically attractive
619: Remarks: computes the gaussian curvatures at each vertex on the boundary of the mesh
622: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
623: */
628: ALE::Obj<ALE::Mesh::sieve_type::supportSet> Hierarchy_curvatureSet(ALE::Obj<ALE::Mesh> original_mesh) {
629: //given a 2D mesh;
630: ALE::Obj<ALE::Mesh::sieve_type::supportSet> output_set;
631: return output_set;
632: }
633: /*=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
634: Hierarchy_defineSpacingFunction
637: Use:
639: Returns the real_section_type associated with "spacing"; which it defines if not already defined
642: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
643: */
651: ALE::Obj<ALE::Mesh::real_section_type> Hierarchy_defineSpacingFunction(ALE::Obj<ALE::Mesh> m, PetscTruth recalculate = PETSC_FALSE) {
652: if (m->hasRealSection("spacing") && !recalculate) {
653: return m->getRealSection("spacing");
654: }
655: ALE::Obj<ALE::Mesh::sieve_type> s = m->getSieve();
656: ALE::Obj<ALE::Mesh::real_section_type> coordinates = m->getRealSection("coordinates");
657: ALE::Obj<ALE::Mesh::label_sequence> vertices = m->depthStratum(0);
658: ALE::Mesh::label_sequence::iterator v_iter = vertices->begin();
659: ALE::Mesh::label_sequence::iterator v_iter_end = vertices->end();
660: ALE::Obj<ALE::Mesh::real_section_type> spacing;
661: int dim = m->getDimension();
662: double *v_coords = new double[dim];
663: if (!m->hasRealSection("spacing")) {
664: spacing = m->getRealSection("spacing");
665: spacing->setFiberDimension(vertices, 1);
666: m->allocate(spacing);
667: } else {
668: spacing = m->getRealSection("spacing");
669: }
670: while (v_iter != v_iter_end) { //standard nearest-neighbor calculation
671: ALE::Obj<ALE::Mesh::sieve_type::coneSet> neighbors = s->cone(s->support(*v_iter));
672: ALE::Mesh::sieve_type::coneSet::iterator n_iter = neighbors->begin();
673: ALE::Mesh::sieve_type::coneSet::iterator n_iter_end = neighbors->end();
674: PetscMemcpy(v_coords, coordinates->restrictPoint(*v_iter), dim*sizeof(double));
675: double min_dist;
676: bool first = true;
677: while (n_iter != n_iter_end) {
678: if (*n_iter != *v_iter) {
679: double dist = 0.;
680: const double * n_coords = coordinates->restrictPoint(*n_iter);
681: for (int i = 0; i < dim; i++) {
682: dist += (n_coords[i] - v_coords[i])*(n_coords[i] - v_coords[i]);
683: }
684: dist = 0.5*sqrt(dist);
685: //PetscPrintf(m->comm(), "%f\n", dist);
686: if (dist < min_dist || first == true) {
687: min_dist = dist;
688: first = false;
689: }
690: }
691: n_iter++;
692: }
693: //PetscPrintf(m->comm(), "%f\n", min_dist);
694: spacing->updatePoint(*v_iter, &min_dist);
695: v_iter++;
696: }
697: delete [] v_coords;
698: return spacing;
699: }
701: /*
702: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
703: Hierarchy_insertVerticesIntoMesh
705: Use: Inserts a bunch of unaffiliated vertices into a boundary mesh for use with
706: the meshbuilder routines
708: Inputs:
710: bound_mesh: the mesh
711: vertex_set: the set of vertices to be inserted
713: Returns: void
715: Remarks:
717: please please please only insert points that are already defined in the coordinate
718: section.
720: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
721: */
726: void Hierarchy_insertVerticesIntoMesh(ALE::Obj<ALE::Mesh> bound_mesh, ALE::Obj<ALE::Mesh::sieve_type::supportSet> vertex_set) {
727: ALE::Mesh::sieve_type::supportSet::iterator vs_iter = vertex_set->begin();
728: ALE::Mesh::sieve_type::supportSet::iterator vs_iter_end = vertex_set->end();
729: ALE::Obj<ALE::Mesh::sieve_type> s = bound_mesh->getSieve();
730: while (vs_iter != vs_iter_end) {
731: if (!bound_mesh->getSieve()->hasPoint(*vs_iter)) {
732: s->addCapPoint(*vs_iter);
733: }
734: vs_iter++;
735: }
736: bound_mesh->stratify();
737: vs_iter = vertex_set->begin();
738: //PetscPrintf(bound_mesh->comm(), "%d is the depth of the inserted points\n", bound_mesh->depth(*vertex_set->begin()));
739: }
741: /*
742: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
743: Hierarchy_coarsenMesh:
745: Use:
747: Returns a coarser version of a 1D, 2D, or 3D mesh
749: Inputs:
751: original_mesh: an effectively 1D simplicial mesh
752: coarsen_factor: the expansion of the spacing function
753: boundary_mesh: a COMPLETE boundary mesh of effective dimension dim-1
755: Output:
757: A coarsened mesh
759: Remarks:
761: Uses routines dependent on having triangle and tetgen installed
763: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
764: */
769: ALE::Obj<ALE::Mesh> Hierarchy_coarsenMesh(ALE::Obj<ALE::Mesh> original_mesh, double coarsen_factor, ALE::Obj<ALE::Mesh> boundary_mesh = PETSC_NULL) {
770: int dim = original_mesh->getDimension();
771: ALE::Obj<ALE::Mesh::real_section_type> spacing = Hierarchy_defineSpacingFunction(original_mesh);
772: //MPI_Comm comm = original_mesh->comm();
773: ALE::Obj<ALE::Mesh::label_sequence> vertices = original_mesh->depthStratum(0);
774: ALE::Mesh::label_sequence::iterator v_iter = vertices->begin();
775: ALE::Mesh::label_sequence::iterator v_iter_end = vertices->end();
776: ALE::Mesh::point_type maxIndex;
777: bool first_maxIndex = true;
778: while (v_iter != v_iter_end) {
779: if (first_maxIndex) {
780: maxIndex = *v_iter;
781: first_maxIndex = false;
782: } else if (maxIndex < *v_iter) {
783: maxIndex = *v_iter;
784: }
785: v_iter++;
786: }
787: if (dim > 3) throw ALE::Exception("Cannot coarsen meshes in more than 3 dimensions");
788: ALE::Obj<ALE::Mesh::sieve_type> s = original_mesh->getSieve();
789: ALE::Obj<ALE::Mesh> output_mesh;
790: if (dim == 1) {
791: ALE::Obj<ALE::Mesh::sieve_type::supportSet> bound_0;
792: if (!boundary_mesh) { //build the 0D boundary
793: bound_0 = Hierarchy_createEffective0DBoundary(original_mesh);
794: } else { //I know this case will never, ever come up, but whatever.
795: ALE::Obj<ALE::Mesh::label_sequence> boundary_vertices = boundary_mesh->depthStratum(0);
796: bound_0 = new ALE::Mesh::sieve_type::supportSet();
797: ALE::Mesh::label_sequence::iterator bv_iter = boundary_vertices->begin();
798: ALE::Mesh::label_sequence::iterator bv_iter_end = boundary_vertices->end();
799: while (bv_iter != bv_iter_end) {
800: bound_0->insert(*bv_iter);
801: bv_iter++;
802: }
803: }
804: }
805: if (dim == 2) {
806: ALE::Obj<ALE::Mesh::sieve_type::supportSet> bound_0;
807: ALE::Obj<ALE::Mesh::sieve_type::supportSet> curvature_set;
808: ALE::Obj<ALE::Mesh> bound_1;
809: ALE::Obj<ALE::Mesh::sieve_type::supportSet> bound_1_vertices = new ALE::Mesh::sieve_type::supportSet();
810: if (!boundary_mesh) {
811: bound_1 = Hierarchy_createEffective1DBoundary(original_mesh, maxIndex);
812: //PetscPrintf(original_mesh->comm(), "%d vertices, %d edges on the boundary\n", bound_1->depthStratum(0)->size(), bound_1->depthStratum(1)->size());
813: }
814: ALE::Obj<ALE::Mesh::label_sequence> bound_1_vert_label = bound_1->depthStratum(0);
815: ALE::Mesh::label_sequence::iterator bv_iter = bound_1_vert_label->begin();
816: ALE::Mesh::label_sequence::iterator bv_iter_end = bound_1_vert_label->end();
817: while (bv_iter != bv_iter_end) {
818: bound_1_vertices->insert(*bv_iter);
819: bv_iter++;
820: }
821: bound_0 = Hierarchy_createEffective0DBoundary(bound_1);
822: //PetscPrintf(comm, "%d vertices in the 0D effective boundary\n", bound_0->size());
823: //coarsen the interior
824: //come up with the interior set:
825: ALE::Obj<ALE::Mesh::sieve_type::supportSet> interior_set = new ALE::Mesh::sieve_type::supportSet();
826: vertices = original_mesh->depthStratum(0);
827: v_iter = vertices->begin();
828: v_iter_end = vertices->end();
829: while (v_iter != v_iter_end) {
830: if (bound_1_vertices->find(*v_iter) == bound_1_vertices->end()) {
831: interior_set->insert(*v_iter);
832: }
833: v_iter++;
834: }
835: //PetscPrintf(original_mesh->comm(), "%d interior vertices\n", interior_set->size());
836: ALE::Obj<ALE::Mesh::sieve_type::supportSet> coarse_interior = Hierarchy_CoarsenVertexSet(original_mesh, spacing, interior_set, bound_1_vertices, coarsen_factor, maxIndex);
837: //coarsen the boundary
838: //PetscPrintf(original_mesh->comm(), "%d vertices in coarsened interior\n", coarse_interior->size());
839: ALE::Obj<ALE::Mesh::sieve_type::supportSet> coarse_bound_set = Hierarchy_CoarsenVertexSet(original_mesh, spacing, bound_1_vertices, bound_0, coarsen_factor, maxIndex); //coarsen the boundary mesh
840: //insert bound_0 into the coarse bound set
841: ALE::Mesh::sieve_type::supportSet::iterator b0_iter = bound_0->begin();
842: ALE::Mesh::sieve_type::supportSet::iterator b0_iter_end = bound_0->end();
843: while (b0_iter != b0_iter_end) {
844: coarse_bound_set->insert(*b0_iter);
845: b0_iter++;
846: }
847: ALE::Obj<ALE::Mesh> coarse_bound = ALE::Surgery_1D_Coarsen_Mesh(bound_1, coarse_bound_set);
848: PetscPrintf(original_mesh->comm(), "%d v, %d e in new coarse boundary\n", coarse_bound->depthStratum(0)->size(), coarse_bound->depthStratum(1)->size());
849: Hierarchy_insertVerticesIntoMesh(coarse_bound, coarse_interior);
850: PetscPrintf(original_mesh->comm(), "%d v, %d e in the thing we feed triangle\n", coarse_bound->depthStratum(0)->size(), coarse_bound->depthStratum(1)->size());
851: //copy over the "marker" label
852: ALE::Obj<ALE::Mesh::label_type> bound_marker_label = coarse_bound->createLabel("marker");
853: ALE::Obj<ALE::Mesh::label_type> marker = original_mesh->getLabel("marker");
854: vertices = coarse_bound->depthStratum(0);
855: v_iter = vertices->begin();
856: v_iter_end = vertices->end();
857: while (v_iter != v_iter_end) {
858: coarse_bound->setValue(bound_marker_label, *v_iter, original_mesh->getValue(marker, *v_iter));
859: v_iter++;
860: }
861: //generate the mesh
862: //this is screwy... we must do this differently
863: coarse_bound->setDimension(1);
864: output_mesh = ALE::Generator<PETSC_MESH_TYPE>::generateMesh(coarse_bound, (original_mesh->depth() != 1), true);
865: output_mesh->stratify();
866: PetscPrintf(original_mesh->comm(), "%d v, %d cells in the output mesh\n", output_mesh->depthStratum(0)->size(), output_mesh->heightStratum(0)->size());
867: }
868: if (dim == 3) {
869: ALE::Obj<ALE::Mesh::sieve_type::supportSet> bound_0;
870: ALE::Obj<ALE::Mesh> bound_1;
871: ALE::Obj<ALE::Mesh::sieve_type::supportSet> bound_1_vertices = new ALE::Mesh::sieve_type::supportSet();
872: ALE::Obj<ALE::Mesh> bound_2;
873: ALE::Obj<ALE::Mesh::sieve_type::supportSet> bound_2_vertices = new ALE::Mesh::sieve_type::supportSet();
874: ALE::Obj<ALE::Mesh::sieve_type::supportSet> interior_vertices = new ALE::Mesh::sieve_type::supportSet();
875: if (!boundary_mesh) {
876: bound_2 = Hierarchy_createEffective2DBoundary(original_mesh, maxIndex);
877: } else {
878: bound_2 = boundary_mesh;
879: }
880: //build the boundary and interior vertex sets from this
881: //PetscPrintf(original_mesh->comm(), "%d vertices, %d cells in 2D boundary\n", bound_2->depthStratum(0)->size(), bound_2->heightStratum(0)->size());
882: ALE::Obj<ALE::Mesh::label_sequence> bound_2_vert_label = bound_2->depthStratum(0);
883: ALE::Mesh::label_sequence::iterator bv_iter = bound_2_vert_label->begin();
884: ALE::Mesh::label_sequence::iterator bv_iter_end = bound_2_vert_label->end();
885: while (bv_iter != bv_iter_end) {
886: bound_2_vertices->insert(*bv_iter);
887: bv_iter++;
888: }
889: v_iter = vertices->begin();
890: v_iter_end = vertices->end();
891: while (v_iter != v_iter_end) {
892: if (bound_2_vertices->find(*v_iter) == bound_2_vertices->end()) {
893: interior_vertices->insert(*v_iter);
894: }
895: v_iter++;
896: }
897: //PetscPrintf(original_mesh->comm(), "%d vertices on the interior\n", interior_vertices->size());
898: bound_1 = Hierarchy_createEffective1DBoundary(bound_2, maxIndex);
899: //PetscPrintf(original_mesh->comm(), "%d vertices, %d edges in 1D boundary\n", bound_1->depthStratum(0)->size(), bound_1->heightStratum(0)->size());
900: ALE::Obj<ALE::Mesh::label_sequence> bound_1_vert_label = bound_1->depthStratum(0);
901: bv_iter = bound_1_vert_label->begin();
902: bv_iter_end = bound_1_vert_label->end();
903: while (bv_iter != bv_iter_end) {
904: bound_1_vertices->insert(*bv_iter);
905: bv_iter++;
906: }
907: bound_0 = Hierarchy_createEffective0DBoundary(bound_1);
908: //PetscPrintf(comm, "%d vertices in the 0D effective boundary\n", bound_0->size());
909: //ok. Coarsen the interior set.
910: ALE::Obj<ALE::Mesh::sieve_type::supportSet> coarse_interior_set = Hierarchy_CoarsenVertexSet(original_mesh, spacing, interior_vertices, bound_2_vertices, coarsen_factor, maxIndex);
911: //PetscPrintf(original_mesh->comm(), "%d vertices in new interior\n", coarse_interior_set->size());
912: //coarsen the boundary mesh
913: ALE::Obj<ALE::Mesh::sieve_type::supportSet> coarse_bound_2_set = Hierarchy_CoarsenVertexSet(original_mesh, spacing, bound_2_vertices, bound_1_vertices, coarsen_factor, maxIndex);
914: //PetscPrintf(original_mesh->comm(), "%d vertices in new boundary interior\n", coarse_bound_2_set->size());
915: //coarsen the boundary skeleton
916: ALE::Obj<ALE::Mesh::sieve_type::supportSet> coarse_bound_1_set = Hierarchy_CoarsenVertexSet(original_mesh, spacing, bound_1_vertices, bound_0, coarsen_factor, maxIndex);
917: //PetscPrintf(original_mesh->comm(), "%d vertices in new boundary skeleton\n", coarse_bound_1_set->size());
918: //merge the coarse_bound_1, bound_0, and coarse_bound_2 sets
919: ALE::Mesh::sieve_type::supportSet::iterator av_iter = coarse_bound_1_set->begin();
920: ALE::Mesh::sieve_type::supportSet::iterator av_iter_end = coarse_bound_1_set->end();
921: while (av_iter != av_iter_end) {
922: coarse_bound_2_set->insert(*av_iter);
923: av_iter++;
924: }
925: av_iter = bound_0->begin();
926: av_iter_end = bound_0->end();
927: while (av_iter != av_iter_end) {
928: coarse_bound_2_set->insert(*av_iter);
929: av_iter++;
930: }
931: //PetscPrintf(original_mesh->comm(), "%d total vertices in the coarse boundary: coarsening by flip.\n", coarse_bound_2_set->size());
932: //ok. Moment of truth; use the FLIPPING to coarsen the boundary mesh
933: ALE::Obj<ALE::Mesh> coarse_bound_mesh = ALE::Surgery_2D_Coarsen_Mesh(bound_2, coarse_bound_2_set, bound_1);
934: //PetscPrintf(original_mesh->comm(), "%d vertices, %d faces in new exterior boundary\n", coarse_bound_mesh->depthStratum(0)->size(), coarse_bound_mesh->heightStratum(0)->size());
935: Hierarchy_insertVerticesIntoMesh(coarse_bound_mesh, coarse_interior_set);
936: //PetscPrintf(original_mesh->comm(), "%d vertices, %d faces to tetgen\n", coarse_bound_mesh->depthStratum(0)->size(), coarse_bound_mesh->heightStratum(0)->size());
937: ALE::Obj<ALE::Mesh::label_type> bound_marker_label = coarse_bound_mesh->createLabel("marker");
938: ALE::Obj<ALE::Mesh::label_type> marker = original_mesh->getLabel("marker");
939: vertices = coarse_bound_mesh->depthStratum(0);
940: v_iter = vertices->begin();
941: v_iter_end = vertices->end();
942: while (v_iter != v_iter_end) {
943: coarse_bound_mesh->setValue(bound_marker_label, *v_iter, original_mesh->getValue(marker, *v_iter));
944: v_iter++;
945: }
946: //generate the mesh
947: //this is screwy... we must do this differently
948: coarse_bound_mesh->setDimension(2);
949: output_mesh = ALE::Generator<PETSC_MESH_TYPE>::generateMesh(coarse_bound_mesh, (original_mesh->depth() != 1), true);
950: output_mesh->stratify();
951: //PetscPrintf(original_mesh->comm(), "%d v, %d cells in the output mesh\n", output_mesh->depthStratum(0)->size(), output_mesh->heightStratum(0)->size());
952: }
953: //PetscPrintf(comm, "leaving overall coarsening\n");
954: return output_mesh;
955: }
957: /*
958: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
959: Hierarchy_createHierarchy:
961: Use:
963: Returns a coarser hierarchy of a 1D, 2D, or 3D mesh
965: Inputs:
967: original_mesh: an effectively 1D simplicial mesh
968: nLevels: the number of levels (including the fine level)
969: coarsen_factor: the expansion of the spacing function
970: boundary_mesh: an optional boundary mesh
971: CtF: do coarse-to-fine hierarchy building (O(ln)) instead of fine-to-coarse (O(n))
973: Output:
975: A coarsened mesh hierarchy
977: Remarks:
979: Uses routines dependent on having triangle and tetgen installed
981: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
982: */
984: ALE::Obj<ALE::Mesh> * Hierarchy_createHierarchy(ALE::Obj<ALE::Mesh> original_mesh,
985: int nLevels,
986: double beta,
987: ALE::Obj<ALE::Mesh> boundary_mesh = PETSC_NULL,
988: PetscTruth CtF = PETSC_FALSE) {
990: ALE::Obj<ALE::Mesh> * meshes = new ALE::Obj<ALE::Mesh>[nLevels];
991: if (CtF) {
992: throw ALE::Exception("Coarse-to-Fine not yet implemented.");
993: } else {
994: //fine-to-coarse hierarchy creation
995: meshes[0] = original_mesh;
996: for (int i = 1; i < nLevels; i++) {
997: meshes[i] = Hierarchy_coarsenMesh(meshes[i-1], beta);
998: }
999: }
1000: return meshes;
1001: }
1003: ALE::Obj<ALE::Mesh> * Hierarchy_createHierarchy_adaptive(ALE::Obj<ALE::Mesh> original_mesh,
1004: unsigned int nVertices,
1005: unsigned int max_meshes,
1006: double beta,
1007: int * nMeshes,
1008: ALE::Obj<ALE::Mesh> boundary_mesh = PETSC_NULL,
1009: PetscTruth CtF = PETSC_FALSE) {
1010: ALE::Obj<ALE::Mesh> * meshes;
1011: if (CtF) {
1012: throw ALE::Exception ("Coarse-to-Fine not yet implemented.");
1013: } else {
1014: //fine to coarse hierarchy creation AS LONG AS the mesh size is over nVertices.
1015: std::list<ALE::Obj<ALE::Mesh> > meshes_list;
1016: meshes_list.clear();
1017: meshes_list.push_front(original_mesh);
1018: ALE::Obj<ALE::Mesh> current_mesh = original_mesh;
1019: while (current_mesh->depthStratum(0)->size() > nVertices && meshes_list.size() <= max_meshes) {
1020: current_mesh = Hierarchy_coarsenMesh(current_mesh, beta);
1021: meshes_list.push_back(current_mesh);
1022: }
1023: meshes = new ALE::Obj<ALE::Mesh>[meshes_list.size()];
1024: *nMeshes = meshes_list.size();
1025: int i = 0;
1026: for (std::list<ALE::Obj<ALE::Mesh> >::iterator m_iter = meshes_list.begin(); m_iter != meshes_list.end(); m_iter++) {
1027: meshes[i] = *m_iter;
1028: i++;
1029: }
1030: }
1031: return meshes;
1032: }
1033: /*
1034: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
1035: Hierarchy_CellsCollide:
1037: Use: tells if the coordinates given by a and b collide
1039: Inputs: a (simplex vertex coords), b (simplex vertex coodrs), dim
1041: Outputs: bool; if they collide or not
1043: Remarks: Sometimes has false positives for shared faces, not a big deal
1045: Assume that
1047: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
1048: */
1050: PetscTruth Hierarchy_CellsCollide(const int dim, const double * a_coords, const double * b_coords) {
1051: int corners = dim+1; //simplices
1052: double * hyperplane = new double[dim];
1053: double * pivotpoint = new double[dim];
1054: double dot_product;
1055: double norm_2;
1056: double dist;
1057: double min_dist = 1.e40;
1058: double max_dist = 0.;
1059: double tolerance = 1.e-30;
1060: bool collides = false;
1061: //initialize the pivot and hyperplane
1062: for (int a = 0; a < corners; a++) {
1063: for (int b = 0; b < corners; b++) {
1064: //find the minimum distance between any two vertices
1065: dist = 0.;
1066: for (int i = 0; i < dim; i++) {
1067: dist += (a_coords[dim*a+i] - b_coords[dim*b+i])*(a_coords[dim*a+i] - b_coords[dim*b+i]);
1068: }
1069: dist = sqrt(dist);
1070: if (dist < min_dist) {
1071: min_dist = dist;
1072: for (int i = 0; i < dim; i++) pivotpoint[i] = 0.5*(a_coords[dim*a+i] + b_coords[dim*b+i]);
1073: }
1074: if (dist > max_dist) {
1075: max_dist = dist;
1076: for (int i = 0; i < dim; i++) hyperplane[i] = (a_coords[dim*a+i] - b_coords[dim*b+i]);
1077: }
1078: }
1079: }
1080: if (max_dist < tolerance) return PETSC_TRUE; //this cell is messed up
1081: //normalize the hyperplane
1082: for (int i = 0; i < dim; i++) {
1083: dist += hyperplane[i]*hyperplane[i];
1084: }
1085: dist = sqrt(dist);
1086: for (int i = 0; i < dim; i++) {
1087: hyperplane[i] = hyperplane[i]/dist;
1088: }
1089: //training run: attempt to properly classify every fine and coarse point; misclassified points get put on the updated hyperplane
1090: //a case; dot products should be positive
1091: for (int a = 0; a < dim+1; a++) {
1092: dot_product = 0.;
1093: norm_2 = 0.;
1094: for (int i = 0; i < dim; i++) {
1095: dot_product += hyperplane[i]*(a_coords[a*dim+i] - pivotpoint[i]);
1096: norm_2 += (a_coords[a*dim+i] - pivotpoint[i])*(a_coords[a*dim+i] - pivotpoint[i]);
1097: }
1098: if (dot_product < 0.-tolerance) { //misclassification
1099: if (fabs(sqrt(norm_2) + dot_product) < tolerance) { //orthogonal direct misclassification case; multiple happenings of this indicate no separator (sorta)
1100: for (int i = 0; i < dim; i++) hyperplane[i] = 0. - hyperplane[i];
1101: } else { //simple misclassification case, just remove the offending component of the hyperplane
1102: for (int i = 0; i < dim; i++) hyperplane[i] -= dot_product*(a_coords[a*dim+i] - pivotpoint[i])/norm_2;
1103: for (int i = 0; i < dim; i++) {
1104: dist += hyperplane[i]*hyperplane[i];
1105: }
1106: dist = sqrt(dist);
1107: for (int i = 0; i < dim; i++) {
1108: hyperplane[i] = hyperplane[i]/dist;
1109: }
1110: }
1111: }
1112: }
1113: //b case; dot products should be negative
1114: for (int b = 0; b < dim+1; b++) {
1115: dot_product = 0.;
1116: norm_2 = 0.;
1117: for (int i = 0; i < dim; i++) {
1118: dot_product += hyperplane[i]*(b_coords[b*dim+i] - pivotpoint[i]);
1119: norm_2 += (b_coords[b*dim+i] - pivotpoint[i])*(b_coords[b*dim+i] - pivotpoint[i]);
1120: }
1121: if (dot_product > tolerance) { //misclassification
1122: if (fabs(sqrt(norm_2) - dot_product) < tolerance) { //orthogonal direct misclassification case; multiple happenings of this indicate no separator (sorta)
1123: for (int i = 0; i < dim; i++) hyperplane[i] = 0. - hyperplane[i];
1124: } else { //simple misclassification case, just remove the offending component of the hyperplane
1125: for (int i = 0; i < dim; i++) hyperplane[i] -= dot_product*(b_coords[b*dim+i] - pivotpoint[i])/norm_2;
1126: for (int i = 0; i < dim; i++) {
1127: dist += hyperplane[i]*hyperplane[i];
1128: }
1129: dist = sqrt(dist);
1130: for (int i = 0; i < dim; i++) {
1131: hyperplane[i] = hyperplane[i]/dist;
1132: }
1133: }
1134: }
1135: }
1136: //testing step: if any vertex is misclassified in this case, we know the figures are intersecting (within tolerance)
1137: for (int a = 0; a < dim+1; a++) {
1138: dot_product = 0.;
1139: for (int i = 0; i < dim; i++) {
1140: dot_product += hyperplane[i]*(a_coords[a*dim+i] - pivotpoint[i]);
1141: }
1142: if (dot_product < 0. - tolerance) collides = true;
1143: }
1144: for (int b = 0; b < dim+1; b++) {
1145: dot_product = 0.;
1146: for (int i = 0; i < dim; i++) {
1147: dot_product += hyperplane[i]*(b_coords[b*dim+i] - pivotpoint[i]);
1148: }
1149: if (dot_product > tolerance) collides = true;
1150: }
1151: delete hyperplane;
1152: delete pivotpoint;
1153: if (collides) return PETSC_TRUE;
1154: return PETSC_FALSE;
1155: }
1156: /*
1157: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
1158: Hierarchy_BBoxesCollide
1160: Sees whether the bounding boxes of the simplices collide; quick check for traversal
1162: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
1163: */
1165: PetscTruth Hierarchy_BBoxesCollide (int dim, const double * a_coords, const double * b_coords) {
1166: double * max_a_coords = new double[dim];
1167: double * min_a_coords = new double[dim];
1168: double * max_b_coords = new double[dim];
1169: double * min_b_coords = new double[dim];
1170: //initialize the maxes and mins
1171: for (int i = 0; i < dim; i++) {
1172: max_a_coords[i] = a_coords[i];
1173: max_b_coords[i] = b_coords[i];
1174: min_a_coords[i] = a_coords[i];
1175: min_b_coords[i] = b_coords[i];
1176: }
1177: //find the min and max in every direction;
1178: for (int i = 0; i < dim; i++) {
1179: for (int v = 0; v < dim+1; v++) {
1180: if (a_coords[v*dim+i] > max_a_coords[i]) max_a_coords[i] = a_coords[v*dim+i];
1181: if (b_coords[v*dim+i] > max_b_coords[i]) max_b_coords[i] = b_coords[v*dim+i];
1182: if (a_coords[v*dim+i] < min_a_coords[i]) min_a_coords[i] = a_coords[v*dim+i];
1183: if (b_coords[v*dim+i] < min_b_coords[i]) min_b_coords[i] = b_coords[v*dim+i];
1184: }
1185: }
1186:
1187: for (int i = 0; i < dim; i++) {
1188: if (max_a_coords[i] < min_a_coords[i]) throw ALE::Exception("a erroneous");
1189: if (max_b_coords[i] < min_b_coords[i]) throw ALE::Exception("b erroneous");
1190: }
1191: PetscTruth collides = PETSC_TRUE;
1192: //if any set of max and min coords are in proper order, then they don't collide;
1193: for (int i = 0; i < dim; i++) {
1194: if ((max_b_coords[i] < min_a_coords[i]) || (max_a_coords[i] < min_b_coords[i])) {
1195: collides = PETSC_FALSE;
1196: }
1197: }
1198: delete max_a_coords;
1199: delete min_a_coords;
1200: delete max_b_coords;
1201: delete min_b_coords;
1202: return collides;
1203: }
1205: /*
1206: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
1207: Hierarchy_qualityInfo:
1209: Use:
1211: Gives quality information on the provided mesh hierarchy
1213: Inputs:
1215: the mesh hierarchy
1216: the number of levels
1218: Output:
1220: prints:
1221: - the min, max, and average of the ratio between circumsphere and insphere radii of tetrahedra per-level
1222: - the ratio between the minimum and maximum radii in the mesh
1223: - the the min, max, and average number of triangles
1224: - the min, max and average ratio of the number of cells between levels
1225: - the min, max, and average of the number of unknowns between levels
1227: =+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=
1228: */
1232: void Hierarchy_qualityInfo(ALE::Obj<ALE::Mesh> * meshes, int nLevels) {
1235: //double min_cells_ratio = 1.e37, max_cells_ratio = 0.;
1236: //double min_unknowns_ratio = 1.e37, max_unknowns_ratio = 0.;
1237: //double tolerance = 1.e-10;
1238: int dim = meshes[0]->getDimension();
1239: double * coords = new double[dim*(dim+1)];
1240: double * fcoords = new double[dim*(dim+1)];
1241: PetscPrintf(meshes[0]->comm(), "_level____|_cells____|_vertices_|_min_len._|_max_len._|_avg_len._|_len_rat._|_min_asp._|_max_asp._|_avg_asp._|_max_ccs._|_avg_ccs._|\n");
1242: for (int current_level = 0; current_level < nLevels; current_level++) {
1243: ALE::Obj<ALE::Mesh> m = meshes[current_level];
1244: ALE::Obj<ALE::Mesh::real_section_type> coordinates = m->getRealSection("coordinates");
1245: ALE::Obj<ALE::Mesh::label_sequence> cells = m->heightStratum(0);
1246: int nCells = cells->size();
1247: int nVertices = m->depthStratum(0)->size();
1248: double min_length_scale = 1.e37, max_length_scale = 0., avg_length_scale = 0.;
1249: double min_aspect_ratio = 1.e37, max_aspect_ratio = 0., avg_aspect_ratio = 0.;
1250: int max_cell_collisions = 0;
1251: double avg_cell_collisions = 0.;
1252: ALE::Mesh::label_sequence::iterator c_iter = cells->begin();
1253: ALE::Mesh::label_sequence::iterator c_iter_end = cells->end();
1254: while (c_iter != c_iter_end) {
1255: //restrict the coordinates of the closure
1256: //Compute the longest edge and perimeter:
1257: PetscMemcpy(coords, m->restrictClosure(coordinates, *c_iter), sizeof(double)*dim*(dim+1));
1258: double max_cell_edge = 0.;
1259: for (int edge = 0; edge < dim+1; edge++) {
1260: //compute the max edge length
1261: double edge_length = 0;
1262: for (int i = 0; i < dim; i++) {
1263: double present_coord = coords[dim*edge+i] - coords[dim*((edge+1)%(dim+1))+i];
1264: edge_length += present_coord*present_coord;
1265: }
1266: edge_length = sqrt(edge_length);
1267: if (edge_length < min_length_scale) min_length_scale = edge_length;
1268: if (edge_length > max_length_scale) max_length_scale = edge_length;
1269: if (edge_length > max_cell_edge) max_cell_edge = edge_length;
1270: }
1271: avg_length_scale += max_cell_edge;
1272: /*Compute the incircle radius:
1273: 2D: r = 2A(abc) / (L(ab) + L(bc) + L(ca))
1274: 3D: r = 3V(abcd) / (A(abc) + A(bcd) + A(acd) + A(cad))
1275: */
1276: double current_inscribed_radius;
1277: double current_aspect_ratio;
1278: if (dim == 2) {
1279: //parallelpiped area
1280: double area = 0.5*fabs((coords[1*2+0]-coords[0*2+0])*(coords[2*2+1]-coords[0*2+1])
1281: - (coords[1*2+1]-coords[0*2+1])*(coords[2*2+0]-coords[0*2+0]));
1282: double perimeter = (sqrt((coords[1*2+0]-coords[0*2+0])*(coords[1*2+0]-coords[0*2+0])
1283: + (coords[1*2+1]-coords[0*2+1])*(coords[1*2+1]-coords[0*2+1])) +
1284: sqrt((coords[1*2+0]-coords[2*2+0])*(coords[1*2+0]-coords[2*2+0])
1285: + (coords[1*2+1]-coords[2*2+1])*(coords[1*2+1]-coords[2*2+1])) +
1286: sqrt((coords[0*2+0]-coords[2*2+0])*(coords[0*2+0]-coords[2*2+0])
1287: + (coords[0*2+1]-coords[2*2+1])*(coords[0*2+1]-coords[2*2+1])));
1288: current_inscribed_radius = 2.*area/perimeter;
1289: //PetscPrintf(m->comm(), "%f inscribed radius %f area %f perimeter\n", current_inscribed_radius, area, perimeter);
1290: } else if (dim == 3) {
1291: double volume = fabs((coords[1*3+0] - coords[0*3+0])*(coords[2*3+1] - coords[0*3+1])*(coords[3*3+2] - coords[0*3+2]) +
1292: (coords[2*3+0] - coords[0*3+0])*(coords[3*3+1] - coords[0*3+1])*(coords[1*3+2] - coords[0*3+2]) +
1293: (coords[3*3+0] - coords[0*3+0])*(coords[1*3+1] - coords[0*3+1])*(coords[2*3+2] - coords[0*3+2]) -
1294: (coords[1*3+0] - coords[0*3+0])*(coords[3*3+1] - coords[0*3+1])*(coords[2*3+2] - coords[0*3+2]) -
1295: (coords[2*3+0] - coords[0*3+0])*(coords[1*3+1] - coords[0*3+1])*(coords[3*3+2] - coords[0*3+2]) -
1296: (coords[3*3+0] - coords[0*3+0])*(coords[2*3+1] - coords[0*3+1])*(coords[1*3+2] - coords[0*3+2]))/6.;
1297: double area = 0.;
1298: //0,1,2
1299: double area_term_1 = (coords[1*3+0] - coords[0*3+0])*(coords[2*3+1] - coords[0*3+1]) - (coords[1*3+1] - coords[0*3+1])*(coords[2*3+0] - coords[0*3+0]);
1300: double area_term_2 = (coords[1*3+1] - coords[0*3+1])*(coords[2*3+2] - coords[0*3+2]) - (coords[1*3+2] - coords[0*3+2])*(coords[2*3+1] - coords[0*3+1]);
1301: double area_term_3 = (coords[1*3+2] - coords[0*3+2])*(coords[2*3+0] - coords[0*3+0]) - (coords[1*3+0] - coords[0*3+0])*(coords[2*3+2] - coords[0*3+2]);
1302: area += sqrt(area_term_1*area_term_1 + area_term_2*area_term_2 + area_term_3*area_term_3);
1303: //0,2,3
1304: area_term_1 = (coords[2*3+0] - coords[0*3+0])*(coords[3*3+1] - coords[0*3+1]) - (coords[2*3+1] - coords[0*3+1])*(coords[3*3+0] - coords[0*3+0]);
1305: area_term_2 = (coords[2*3+1] - coords[0*3+1])*(coords[3*3+2] - coords[0*3+2]) - (coords[2*3+2] - coords[0*3+2])*(coords[3*3+1] - coords[0*3+1]);
1306: area_term_3 = (coords[2*3+2] - coords[0*3+2])*(coords[3*3+0] - coords[0*3+0]) - (coords[2*3+0] - coords[0*3+0])*(coords[3*3+2] - coords[0*3+2]);
1307: area += 0.5*sqrt(area_term_1*area_term_1 + area_term_2*area_term_2 + area_term_3*area_term_3);
1308: //0,1,3
1309: area_term_1 = (coords[1*3+0] - coords[0*3+0])*(coords[3*3+1] - coords[0*3+1]) - (coords[1*3+1] - coords[0*3+1])*(coords[3*3+0] - coords[0*3+0]);
1310: area_term_2 = (coords[1*3+1] - coords[0*3+1])*(coords[3*3+2] - coords[0*3+2]) - (coords[1*3+2] - coords[0*3+2])*(coords[3*3+1] - coords[0*3+1]);
1311: area_term_3 = (coords[1*3+2] - coords[0*3+2])*(coords[3*3+0] - coords[0*3+0]) - (coords[1*3+0] - coords[0*3+0])*(coords[3*3+2] - coords[0*3+2]);
1312: area += 0.5*sqrt(area_term_1*area_term_1 + area_term_2*area_term_2 + area_term_3*area_term_3);
1313: //1,2,3
1314: area_term_1 = (coords[2*3+0] - coords[1*3+0])*(coords[3*3+1] - coords[1*3+1]) - (coords[2*3+1] - coords[1*3+1])*(coords[3*3+0] - coords[1*3+0]);
1315: area_term_2 = (coords[2*3+1] - coords[1*3+1])*(coords[3*3+2] - coords[1*3+2]) - (coords[2*3+2] - coords[1*3+2])*(coords[3*3+1] - coords[1*3+1]);
1316: area_term_3 = (coords[2*3+2] - coords[1*3+2])*(coords[3*3+0] - coords[1*3+0]) - (coords[2*3+0] - coords[1*3+0])*(coords[3*3+2] - coords[1*3+2]);
1317: area += 0.5*sqrt(area_term_1*area_term_1 + area_term_2*area_term_2 + area_term_3*area_term_3);
1318: current_inscribed_radius = 3.*volume/area;
1319: //PetscPrintf(m->comm(), "%f inscribed radius, %f volume, %f surface area\n", current_inscribed_radius, volume, area);
1320: }
1321: current_aspect_ratio = max_cell_edge/(current_inscribed_radius*2);
1322: if (current_aspect_ratio > max_aspect_ratio) max_aspect_ratio = current_aspect_ratio;
1323: if (current_aspect_ratio < min_aspect_ratio) min_aspect_ratio = current_aspect_ratio;
1324: avg_aspect_ratio += current_aspect_ratio;
1325: //PetscPrintf(m->comm(), "%f aspect ratio\n", current_aspect_ratio);
1326: c_iter++;
1327: } //end of loop over cells
1328: //ok we need to do the standard location thingamajob;
1329: c_iter = cells->begin();
1330: c_iter_end = cells->end();
1331: if (current_level != 0) {
1332: /*
1333: while (c_iter != c_iter_end) {
1334: int cell_collisions = 0;
1335: //now we determine the number of collisions between each tri/tet in a coarse mesh with those in the next finer mesh
1336: //IF THERE IS NO LINEAR SEPARATOR BETWEEN THE TWO, THEN THEY COLLIDE
1337: ALE::Obj<ALE::Mesh> f_m = meshes[current_level-1]; //get the fine mesh
1338: ALE::Obj<ALE::Mesh::label_sequence> f_cells = f_m->heightStratum(0);
1339: ALE::Obj<ALE::Mesh::real_section_type> f_coordinates = f_m->getRealSection("coordinates");
1340: ALE::Mesh::label_sequence::iterator fc_iter = f_cells->begin();
1341: ALE::Mesh::label_sequence::iterator fc_iter_end = f_cells->end();
1342: while (fc_iter != fc_iter_end) {
1343: PetscMemcpy(fcoords, f_m->restrictClosure(f_coordinates, *fc_iter), sizeof(double)*dim*(dim+1));
1344: if (Hierarchy_CellsCollide(dim, coords, fcoords) == PETSC_TRUE) cell_collisions++;
1345: fc_iter++;
1346: }
1347: //PetscPrintf(m->comm(), "cell collisions: %d\n", cell_collisions);
1348: c_iter++;
1349: }
1350: */
1351: //NEW STUFF
1353: ALE::Obj<ALE::Mesh> f_m = meshes[current_level-1]; //get the fine mesh
1354: ALE::Obj<ALE::Mesh::sieve_type> f_s = f_m->getSieve();
1355: ALE::Obj<ALE::Mesh::sieve_type> s = m->getSieve();
1356: ALE::Obj<ALE::Mesh::real_section_type> f_coordinates = f_m->getRealSection("coordinates");
1357: ALE::Obj<ALE::Mesh::label_sequence> fcells = f_m->heightStratum(0);
1359: ALE::Obj<ALE::Mesh::sieve_type::supportSet> c_traversal = new ALE::Mesh::sieve_type::supportSet();
1360: ALE::Obj<ALE::Mesh::sieve_type::supportSet> f_traversal = new ALE::Mesh::sieve_type::supportSet();
1362: std::list<ALE::Mesh::point_type> c_cell_list;
1363: std::list<ALE::Mesh::point_type> f_cell_guesses;
1364: std::list<ALE::Mesh::point_type> f_cell_list;
1366: c_iter = cells->begin();
1367: c_iter_end = cells->end();
1368: c_traversal->clear();
1369: f_traversal->clear();
1370: while(c_iter != c_iter_end) {
1371: PetscMemcpy(coords, m->restrictClosure(coordinates, *c_iter), sizeof(double)*dim*(dim+1));
1372: if (c_traversal->find(*c_iter) == c_traversal->end()) {
1373: //locate an initial colliding cell
1374: ALE::Mesh::label_sequence::iterator f_iter = fcells->begin();
1375: ALE::Mesh::label_sequence::iterator f_iter_end = fcells->end();
1376: bool outer_located = false;
1377: while (f_iter != f_iter_end && !outer_located) {
1378: PetscMemcpy(fcoords, f_m->restrictClosure(f_coordinates, *f_iter), sizeof(double)*dim*(dim+1));
1379: if (Hierarchy_BBoxesCollide(dim, coords, fcoords)) {
1380: outer_located = true;
1381: c_cell_list.push_front(*c_iter);
1382: f_cell_guesses.push_front(*f_iter);
1383: c_traversal->insert(*c_iter);
1384: f_traversal->insert(*f_iter);
1385: }
1386: f_iter++;
1387: }
1388: while (!c_cell_list.empty()) {
1389: int nCollisions = 0;
1390: int nComparisons = 0;
1391: int nBBoxCollisions = 0;
1392: ALE::Mesh::point_type c_current_cell = c_cell_list.front();
1393: c_cell_list.pop_front();
1394: PetscMemcpy(coords, m->restrictClosure(coordinates, c_current_cell), sizeof(double)*dim*(dim+1));
1395: ALE::Mesh::point_type f_current_guess = f_cell_guesses.front();
1396: f_cell_guesses.pop_front();
1397: bool found_the_boundbox = false;
1398: //traverse outward from the fine guess
1399: f_cell_list.push_front(f_current_guess);
1400: f_traversal->insert(f_current_guess);
1401: while (!f_cell_list.empty()) {
1402: nComparisons++;
1403: ALE::Mesh::point_type f_current_cell = f_cell_list.front();
1404: f_cell_list.pop_front();
1405: PetscMemcpy(fcoords, f_m->restrictClosure(f_coordinates, f_current_cell), sizeof(double)*dim*(dim+1));
1406: bool bbox_collide = (Hierarchy_BBoxesCollide(dim, coords, fcoords));
1407: //if we have yet to find the box, then we have an unrestricted search; if we have found the box, then we only search within the box
1408: if (bbox_collide || !found_the_boundbox) {
1409: if (bbox_collide) {
1410: if (!found_the_boundbox) {
1411: //clear the current queue and traversal list
1412: f_traversal->clear();
1413: f_cell_list.clear();
1414: }
1415: found_the_boundbox = true;
1416: f_current_guess = f_current_cell;
1417: nBBoxCollisions++;
1418: }
1419: //test to see if this fine cell ACTUALLY collides with the coarse one
1420: if (bbox_collide) {
1421: if (Hierarchy_CellsCollide(dim, coords, fcoords)) nCollisions++;
1422: }
1423: ALE::Obj<ALE::Mesh::sieve_type::coneSet> fine_neighbors = f_s->support(f_s->cone(f_current_cell));
1424: ALE::Mesh::sieve_type::coneSet::iterator fn_iter = fine_neighbors->begin();
1425: ALE::Mesh::sieve_type::coneSet::iterator fn_iter_end = fine_neighbors->end();
1426: //add the neighbors
1427: while (fn_iter != fn_iter_end) {
1428: if (f_traversal->find(*fn_iter) == f_traversal->end()) {
1429: f_cell_list.push_back(*fn_iter);
1430: f_traversal->insert(*fn_iter);
1431: }
1432: fn_iter++;
1433: }
1434: }
1435: }
1436: //clear the fine traversal list
1437: f_traversal->clear();
1438: //add the neighbors of the coarse cell to the list if they haven't yet been covered
1439: ALE::Obj<ALE::Mesh::sieve_type::coneSet> neighbors = s->support(s->cone(c_current_cell));
1440: ALE::Mesh::sieve_type::coneSet::iterator n_iter = neighbors->begin();
1441: ALE::Mesh::sieve_type::coneSet::iterator n_iter_end = neighbors->end();
1442: while (n_iter != n_iter_end) {
1443: if (c_traversal->find(*n_iter) == c_traversal->end()) {
1444: c_cell_list.push_back(*n_iter);
1445: f_cell_guesses.push_back(f_current_guess);
1446: c_traversal->insert(*n_iter);
1447: }
1448: n_iter++;
1449: }
1450: if (max_cell_collisions < nCollisions) max_cell_collisions = nCollisions;
1451: avg_cell_collisions += nCollisions;
1452: // PetscPrintf(m->comm(), "collisions: %d, bound-box collisions: %d, comparisons %d\n", nCollisions, nBBoxCollisions, nComparisons);
1453: }
1454: }
1455: c_iter++;
1456: }
1457: }
1458: //output the stats for this level:
1459: //PetscPrintf(m->comm(), "_level____|_cells____|_vertices_|_min_len._|_max_len._|_avg_len._|_len_rat._|_min_asp._|_max_asp._|_avg_asp._|_max_ccs._|_avg_ccs._|\n");
1460: avg_length_scale = avg_length_scale / nCells;
1461: avg_aspect_ratio = avg_aspect_ratio / nCells;
1462: avg_cell_collisions = avg_cell_collisions / nCells;
1463: PetscPrintf(m->comm(), " %9d %9d %9d %9f %9f %9f %9f %9f %9f %9f %9d %9f\n", current_level, nCells, nVertices, min_length_scale, max_length_scale, avg_length_scale, min_length_scale / max_length_scale, min_aspect_ratio, max_aspect_ratio, avg_aspect_ratio, max_cell_collisions, avg_cell_collisions);
1465: }
1466: delete coords;
1467: delete fcoords;
1468: }