MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_VisualizationHelpers_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 
47 #ifndef MUELU_VISUALIZATIONHELPERS_DEF_HPP_
48 #define MUELU_VISUALIZATIONHELPERS_DEF_HPP_
49 
50 #include <Xpetra_Matrix.hpp>
51 #include <Xpetra_CrsMatrixWrap.hpp>
52 #include <Xpetra_ImportFactory.hpp>
55 #include "MueLu_Level.hpp"
56 #include "MueLu_Graph.hpp"
57 #include "MueLu_Monitor.hpp"
58 #include <vector>
59 #include <list>
60 #include <algorithm>
61 #include <string>
62 
63 #ifdef HAVE_MUELU_CGAL
64 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
65 #include <CGAL/convex_hull_2.h>
66 #include <CGAL/Polyhedron_3.h>
67 #include <CGAL/convex_hull_3.h>
68 #endif
69 
70 
71 namespace MueLu {
72 
73  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75  RCP<ParameterList> validParamList = rcp(new ParameterList());
76 
77  validParamList->set< std::string > ("visualization: output filename", "viz%LEVELID", "filename for VTK-style visualization output");
78  validParamList->set< int > ("visualization: output file: time step", 0, "time step variable for output file name");// Remove me?
79  validParamList->set< int > ("visualization: output file: iter", 0, "nonlinear iteration variable for output file name");//Remove me?
80  validParamList->set<std::string> ("visualization: style", "Point Cloud", "style of aggregate visualization for VTK output. Can be 'Point Cloud', 'Jacks', 'Convex Hulls'");
81  validParamList->set<bool> ("visualization: build colormap", false, "Whether to build a random color map in a separate xml file.");
82  validParamList->set<bool> ("visualization: fine graph edges", false, "Whether to draw all fine node connections along with the aggregates.");
83 
84  return validParamList;
85  }
86 
87  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
88  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doPointCloud(std::vector<int>& vertices, std::vector<int>& geomSizes, LO /* numLocalAggs */, LO numFineNodes) {
89  vertices.reserve(numFineNodes);
90  geomSizes.reserve(numFineNodes);
91  for(LocalOrdinal i = 0; i < numFineNodes; i++)
92  {
93  vertices.push_back(i);
94  geomSizes.push_back(1);
95  }
96  }
97 
98  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
99  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doJacks(std::vector<int>& vertices, std::vector<int>& geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector<bool>& isRoot, const ArrayRCP<LO>& vertex2AggId) {
100  //For each aggregate, find the root node then connect it with the other nodes in the aggregate
101  //Doesn't matter the order, as long as all edges are found.
102  vertices.reserve(vertices.size() + 3 * (numFineNodes - numLocalAggs));
103  geomSizes.reserve(vertices.size() + 2 * (numFineNodes - numLocalAggs));
104  int root = 0;
105  for(int i = 0; i < numLocalAggs; i++) //TODO: Replace this O(n^2) with a better way
106  {
107  while(!isRoot[root])
108  root++;
109  int numInAggFound = 0;
110  for(int j = 0; j < numFineNodes; j++)
111  {
112  if(j == root) //don't make a connection from the root to itself
113  {
114  numInAggFound++;
115  continue;
116  }
117  if(vertex2AggId[root] == vertex2AggId[j])
118  {
119  vertices.push_back(root);
120  vertices.push_back(j);
121  geomSizes.push_back(2);
122  //Also draw the free endpoint explicitly for the current line
123  vertices.push_back(j);
124  geomSizes.push_back(1);
125  numInAggFound++;
126  //if(numInAggFound == aggSizes_[vertex2AggId_[root]]) //don't spend more time looking if done with that root
127  // break;
128  }
129  }
130  root++; //get set up to look for the next root
131  }
132  }
133 
134  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
135  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doConvexHulls2D(std::vector<int>& vertices, std::vector<int>& geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector<bool>& /* isRoot */, const ArrayRCP<LO>& vertex2AggId, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& xCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& yCoords) {
136  for(int agg = 0; agg < numLocalAggs; agg++) {
137  std::list<int> aggNodes;
138  for(int i = 0; i < numFineNodes; i++) {
139  if(vertex2AggId[i] == agg)
140  aggNodes.push_back(i);
141  }
142  //have a list of nodes in the aggregate
144  "CoarseningVisualization::doConvexHulls2D: aggregate contains zero nodes!");
145  if(aggNodes.size() == 1) {
146  vertices.push_back(aggNodes.front());
147  geomSizes.push_back(1);
148  continue;
149  }
150  if(aggNodes.size() == 2) {
151  vertices.push_back(aggNodes.front());
152  vertices.push_back(aggNodes.back());
153  geomSizes.push_back(2);
154  continue;
155  }
156  //check if all points are collinear, need to explicitly draw a line in that case.
157  bool collinear = true; //assume true at first, if a segment not parallel to others then clear
158  {
159  std::list<int>::iterator it = aggNodes.begin();
160  myVec3 firstPoint(xCoords[*it], yCoords[*it], 0);
161  it++;
162  myVec3 secondPoint(xCoords[*it], yCoords[*it], 0);
163  it++; //it now points to third node in the aggregate
164  myVec3 norm1(-(secondPoint.y - firstPoint.y), secondPoint.x - firstPoint.x, 0);
165  do {
166  myVec3 thisNorm(yCoords[*it] - firstPoint.y, firstPoint.x - xCoords[*it], 0);
167  //rotate one of the vectors by 90 degrees so that dot product is 0 if the two are parallel
168  double temp = thisNorm.x;
169  thisNorm.x = thisNorm.y;
170  thisNorm.y = temp;
171  double comp = dotProduct(norm1, thisNorm);
172  if(-1e-8 > comp || comp > 1e-8) {
173  collinear = false;
174  break;
175  }
176  it++;
177  }
178  while(it != aggNodes.end());
179  }
180  if(collinear)
181  {
182  //find the most distant two points in the plane and use as endpoints of line representing agg
183  std::list<int>::iterator min = aggNodes.begin(); //min X then min Y where x is a tie
184  std::list<int>::iterator max = aggNodes.begin(); //max X then max Y where x is a tie
185  for(std::list<int>::iterator it = ++aggNodes.begin(); it != aggNodes.end(); it++) {
186  if(xCoords[*it] < xCoords[*min])
187  min = it;
188  else if(xCoords[*it] == xCoords[*min]) {
189  if(yCoords[*it] < yCoords[*min])
190  min = it;
191  }
192  if(xCoords[*it] > xCoords[*max])
193  max = it;
194  else if(xCoords[*it] == xCoords[*max]) {
195  if(yCoords[*it] > yCoords[*max])
196  max = it;
197  }
198  }
199  //Just set up a line between nodes *min and *max
200  vertices.push_back(*min);
201  vertices.push_back(*max);
202  geomSizes.push_back(2);
203  continue; //jump to next aggregate in loop
204  }
205  std::vector<myVec2> points;
206  std::vector<int> nodes;
207  for(std::list<int>::iterator it = aggNodes.begin(); it != aggNodes.end(); it++) {
208  points.push_back(myVec2(xCoords[*it], yCoords[*it]));
209  nodes.push_back(*it);
210  }
211  std::vector<int> hull = giftWrap(points, nodes, xCoords, yCoords);
212  vertices.reserve(vertices.size() + hull.size());
213  for(size_t i = 0; i < hull.size(); i++) {
214  vertices.push_back(hull[i]);
215  }
216  geomSizes.push_back(hull.size());
217  }
218  }
219 
220 #ifdef HAVE_MUELU_CGAL
221  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
222  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doCGALConvexHulls2D(std::vector<int>& vertices, std::vector<int>& geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector<bool>& isRoot, const ArrayRCP<LO>& vertex2AggId, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& xCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& yCoords) {
223 
224  typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
225  typedef K::Point_2 Point_2;
226 
227  for(int agg = 0; agg < numLocalAggs; agg++) {
228  std::vector<int> aggNodes;
229  std::vector<Point_2> aggPoints;
230  for(int i = 0; i < numFineNodes; i++) {
231  if(vertex2AggId[i] == agg) {
232  Point_2 p(xCoords[i], yCoords[i]);
233  aggPoints.push_back(p);
234  aggNodes.push_back(i);
235  }
236  }
237  //have a list of nodes in the aggregate
239  "CoarseningVisualization::doCGALConvexHulls2D: aggregate contains zero nodes!");
240  if(aggNodes.size() == 1) {
241  vertices.push_back(aggNodes.front());
242  geomSizes.push_back(1);
243  continue;
244  }
245  if(aggNodes.size() == 2) {
246  vertices.push_back(aggNodes.front());
247  vertices.push_back(aggNodes.back());
248  geomSizes.push_back(2);
249  continue;
250  }
251  //check if all points are collinear, need to explicitly draw a line in that case.
252  bool collinear = true; //assume true at first, if a segment not parallel to others then clear
253  {
254  std::vector<int>::iterator it = aggNodes.begin();
255  myVec3 firstPoint(xCoords[*it], yCoords[*it], 0);
256  it++;
257  myVec3 secondPoint(xCoords[*it], yCoords[*it], 0);
258  it++; //it now points to third node in the aggregate
259  myVec3 norm1(-(secondPoint.y - firstPoint.y), secondPoint.x - firstPoint.x, 0);
260  do {
261  myVec3 thisNorm(yCoords[*it] - firstPoint.y, firstPoint.x - xCoords[*it], 0);
262  //rotate one of the vectors by 90 degrees so that dot product is 0 if the two are parallel
263  double temp = thisNorm.x;
264  thisNorm.x = thisNorm.y;
265  thisNorm.y = temp;
266  double comp = dotProduct(norm1, thisNorm);
267  if(-1e-8 > comp || comp > 1e-8) {
268  collinear = false;
269  break;
270  }
271  it++;
272  }
273  while(it != aggNodes.end());
274  }
275  if(collinear)
276  {
277  //find the most distant two points in the plane and use as endpoints of line representing agg
278  std::vector<int>::iterator min = aggNodes.begin(); //min X then min Y where x is a tie
279  std::vector<int>::iterator max = aggNodes.begin(); //max X then max Y where x is a tie
280  for(std::vector<int>::iterator it = ++aggNodes.begin(); it != aggNodes.end(); it++) {
281  if(xCoords[*it] < xCoords[*min])
282  min = it;
283  else if(xCoords[*it] == xCoords[*min]) {
284  if(yCoords[*it] < yCoords[*min])
285  min = it;
286  }
287  if(xCoords[*it] > xCoords[*max])
288  max = it;
289  else if(xCoords[*it] == xCoords[*max]) {
290  if(yCoords[*it] > yCoords[*max])
291  max = it;
292  }
293  }
294  //Just set up a line between nodes *min and *max
295  vertices.push_back(*min);
296  vertices.push_back(*max);
297  geomSizes.push_back(2);
298  continue; //jump to next aggregate in loop
299  }
300  // aggregate has more than 2 points and is not collinear
301  {
302  std::vector<Point_2> result;
303  CGAL::convex_hull_2( aggPoints.begin(), aggPoints.end(), std::back_inserter(result) );
304  // loop over all result points and find the corresponding node id
305  for (size_t r = 0; r < result.size(); r++) {
306  // loop over all aggregate nodes and find corresponding node id
307  for(size_t l = 0; l < aggNodes.size(); l++)
308  {
309  if(fabs(result[r].x() - xCoords[aggNodes[l]]) < 1e-12 &&
310  fabs(result[r].y() - yCoords[aggNodes[l]]) < 1e-12)
311  {
312  vertices.push_back(aggNodes[l]);
313  break;
314  }
315  }
316 
317  }
318  geomSizes.push_back(result.size());
319  }
320  }
321  }
322 
323 #endif
324 
325  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doConvexHulls3D(std::vector<int>& vertices, std::vector<int>& geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector<bool>& /* isRoot */, const ArrayRCP<LO>& vertex2AggId, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& xCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& yCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& zCoords) {
327  //Use 3D quickhull algo.
328  //Vector of node indices representing triangle vertices
329  //Note: Calculate the hulls first since will only include point data for points in the hulls
330  //Effectively the size() of vertIndices after each hull is added to it
331  typedef std::list<int>::iterator Iter;
332  for(int agg = 0; agg < numLocalAggs; agg++) {
333  std::list<int> aggNodes; //At first, list of all nodes in the aggregate. As nodes are enclosed or included by/in hull, remove them
334  for(int i = 0; i < numFineNodes; i++) {
335  if(vertex2AggId[i] == agg)
336  aggNodes.push_back(i);
337  }
338  //First, check anomalous cases
340  "CoarseningVisualization::doConvexHulls3D: aggregate contains zero nodes!");
341  if(aggNodes.size() == 1) {
342  vertices.push_back(aggNodes.front());
343  geomSizes.push_back(1);
344  continue;
345  } else if(aggNodes.size() == 2) {
346  vertices.push_back(aggNodes.front());
347  vertices.push_back(aggNodes.back());
348  geomSizes.push_back(2);
349  continue;
350  }
351  //check for collinearity
352  bool areCollinear = true;
353  {
354  Iter it = aggNodes.begin();
355  myVec3 firstVec(xCoords[*it], yCoords[*it], zCoords[*it]);
356  myVec3 comp;
357  {
358  it++;
359  myVec3 secondVec(xCoords[*it], yCoords[*it], zCoords[*it]); //cross this with other vectors to compare
360  comp = vecSubtract(secondVec, firstVec);
361  it++;
362  }
363  for(; it != aggNodes.end(); it++) {
364  myVec3 thisVec(xCoords[*it], yCoords[*it], zCoords[*it]);
365  myVec3 cross = crossProduct(vecSubtract(thisVec, firstVec), comp);
366  if(mymagnitude(cross) > 1e-10) {
367  areCollinear = false;
368  break;
369  }
370  }
371  }
372  if(areCollinear)
373  {
374  //find the endpoints of segment describing all the points
375  //compare x, if tie compare y, if tie compare z
376  Iter min = aggNodes.begin();
377  Iter max = aggNodes.begin();
378  Iter it = ++aggNodes.begin();
379  for(; it != aggNodes.end(); it++) {
380  if(xCoords[*it] < xCoords[*min]) min = it;
381  else if(xCoords[*it] == xCoords[*min]) {
382  if(yCoords[*it] < yCoords[*min]) min = it;
383  else if(yCoords[*it] == yCoords[*min]) {
384  if(zCoords[*it] < zCoords[*min]) min = it;
385  }
386  }
387  if(xCoords[*it] > xCoords[*max]) max = it;
388  else if(xCoords[*it] == xCoords[*max]) {
389  if(yCoords[*it] > yCoords[*max]) max = it;
390  else if(yCoords[*it] == yCoords[*max]) {
391  if(zCoords[*it] > zCoords[*max])
392  max = it;
393  }
394  }
395  }
396  vertices.push_back(*min);
397  vertices.push_back(*max);
398  geomSizes.push_back(2);
399  continue;
400  }
401  bool areCoplanar = true;
402  {
403  //number of points is known to be >= 3
404  Iter vert = aggNodes.begin();
405  myVec3 v1(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
406  vert++;
407  myVec3 v2(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
408  vert++;
409  myVec3 v3(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
410  vert++;
411  //Make sure the first three points aren't also collinear (need a non-degenerate triangle to get a normal)
412  while(mymagnitude(crossProduct(vecSubtract(v1, v2), vecSubtract(v1, v3))) < 1e-10) {
413  //Replace the third point with the next point
414  v3 = myVec3(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
415  vert++;
416  }
417  for(; vert != aggNodes.end(); vert++) {
418  myVec3 pt(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
419  if(fabs(pointDistFromTri(pt, v1, v2, v3)) > 1e-12) {
420  areCoplanar = false;
421  break;
422  }
423  }
424  if(areCoplanar) {
425  //do 2D convex hull
426  myVec3 planeNorm = getNorm(v1, v2, v3);
427  planeNorm.x = fabs(planeNorm.x);
428  planeNorm.y = fabs(planeNorm.y);
429  planeNorm.z = fabs(planeNorm.z);
430  std::vector<myVec2> points;
431  std::vector<int> nodes;
432  if(planeNorm.x >= planeNorm.y && planeNorm.x >= planeNorm.z) {
433  //project points to yz plane to make hull
434  for(Iter it = aggNodes.begin(); it != aggNodes.end(); it++) {
435  nodes.push_back(*it);
436  points.push_back(myVec2(yCoords[*it], zCoords[*it]));
437  }
438  }
439  if(planeNorm.y >= planeNorm.x && planeNorm.y >= planeNorm.z) {
440  //use xz
441  for(Iter it = aggNodes.begin(); it != aggNodes.end(); it++) {
442  nodes.push_back(*it);
443  points.push_back(myVec2(xCoords[*it], zCoords[*it]));
444  }
445  }
446  if(planeNorm.z >= planeNorm.x && planeNorm.z >= planeNorm.y) {
447  for(Iter it = aggNodes.begin(); it != aggNodes.end(); it++) {
448  nodes.push_back(*it);
449  points.push_back(myVec2(xCoords[*it], yCoords[*it]));
450  }
451  }
452  std::vector<int> convhull2d = giftWrap(points, nodes, xCoords, yCoords);
453  geomSizes.push_back(convhull2d.size());
454  vertices.reserve(vertices.size() + convhull2d.size());
455  for(size_t i = 0; i < convhull2d.size(); i++)
456  vertices.push_back(convhull2d[i]);
457  continue;
458  }
459  }
460  Iter exIt = aggNodes.begin(); //iterator to be used for searching for min/max x/y/z
461  int extremeSix[] = {*exIt, *exIt, *exIt, *exIt, *exIt, *exIt}; //nodes with minimumX, maxX, minY ...
462  exIt++;
463  for(; exIt != aggNodes.end(); exIt++) {
464  Iter it = exIt;
465  if(xCoords[*it] < xCoords[extremeSix[0]] ||
466  (xCoords[*it] == xCoords[extremeSix[0]] && yCoords[*it] < yCoords[extremeSix[0]]) ||
467  (xCoords[*it] == xCoords[extremeSix[0]] && yCoords[*it] == yCoords[extremeSix[0]] && zCoords[*it] < zCoords[extremeSix[0]]))
468  extremeSix[0] = *it;
469  if(xCoords[*it] > xCoords[extremeSix[1]] ||
470  (xCoords[*it] == xCoords[extremeSix[1]] && yCoords[*it] > yCoords[extremeSix[1]]) ||
471  (xCoords[*it] == xCoords[extremeSix[1]] && yCoords[*it] == yCoords[extremeSix[1]] && zCoords[*it] > zCoords[extremeSix[1]]))
472  extremeSix[1] = *it;
473  if(yCoords[*it] < yCoords[extremeSix[2]] ||
474  (yCoords[*it] == yCoords[extremeSix[2]] && zCoords[*it] < zCoords[extremeSix[2]]) ||
475  (yCoords[*it] == yCoords[extremeSix[2]] && zCoords[*it] == zCoords[extremeSix[2]] && xCoords[*it] < xCoords[extremeSix[2]]))
476  extremeSix[2] = *it;
477  if(yCoords[*it] > yCoords[extremeSix[3]] ||
478  (yCoords[*it] == yCoords[extremeSix[3]] && zCoords[*it] > zCoords[extremeSix[3]]) ||
479  (yCoords[*it] == yCoords[extremeSix[3]] && zCoords[*it] == zCoords[extremeSix[3]] && xCoords[*it] > xCoords[extremeSix[3]]))
480  extremeSix[3] = *it;
481  if(zCoords[*it] < zCoords[extremeSix[4]] ||
482  (zCoords[*it] == zCoords[extremeSix[4]] && xCoords[*it] < xCoords[extremeSix[4]]) ||
483  (zCoords[*it] == zCoords[extremeSix[4]] && xCoords[*it] == xCoords[extremeSix[4]] && yCoords[*it] < yCoords[extremeSix[4]]))
484  extremeSix[4] = *it;
485  if(zCoords[*it] > zCoords[extremeSix[5]] ||
486  (zCoords[*it] == zCoords[extremeSix[5]] && xCoords[*it] > xCoords[extremeSix[5]]) ||
487  (zCoords[*it] == zCoords[extremeSix[5]] && xCoords[*it] == xCoords[extremeSix[5]] && yCoords[*it] > zCoords[extremeSix[5]]))
488  extremeSix[5] = *it;
489  }
490  myVec3 extremeVectors[6];
491  for(int i = 0; i < 6; i++) {
492  myVec3 thisExtremeVec(xCoords[extremeSix[i]], yCoords[extremeSix[i]], zCoords[extremeSix[i]]);
493  extremeVectors[i] = thisExtremeVec;
494  }
495  double maxDist = 0;
496  int base1 = 0; //ints from 0-5: which pair out of the 6 extreme points are the most distant? (indices in extremeSix and extremeVectors)
497  int base2 = 0;
498  for(int i = 0; i < 5; i++) {
499  for(int j = i + 1; j < 6; j++) {
500  double thisDist = distance(extremeVectors[i], extremeVectors[j]);
501  if(thisDist > maxDist) {
502  maxDist = thisDist;
503  base1 = i;
504  base2 = j;
505  }
506  }
507  }
508  std::list<myTriangle> hullBuilding; //each Triangle is a triplet of nodes (int IDs) that form a triangle
509  //remove base1 and base2 iters from aggNodes, they are known to be in the hull
510  aggNodes.remove(extremeSix[base1]);
511  aggNodes.remove(extremeSix[base2]);
512  //extremeSix[base1] and [base2] still have the myVec3 representation
513  myTriangle tri;
514  tri.v1 = extremeSix[base1];
515  tri.v2 = extremeSix[base2];
516  //Now find the point that is furthest away from the line between base1 and base2
517  maxDist = 0;
518  //need the vectors to do "quadruple product" formula
519  myVec3 b1 = extremeVectors[base1];
520  myVec3 b2 = extremeVectors[base2];
521  Iter thirdNode;
522  for(Iter node = aggNodes.begin(); node != aggNodes.end(); node++) {
523  myVec3 nodePos(xCoords[*node], yCoords[*node], zCoords[*node]);
524  double dist = mymagnitude(crossProduct(vecSubtract(nodePos, b1), vecSubtract(nodePos, b2))) / mymagnitude(vecSubtract(b2, b1));
525  if(dist > maxDist) {
526  maxDist = dist;
527  thirdNode = node;
528  }
529  }
530  //Now know the last node in the first triangle
531  tri.v3 = *thirdNode;
532  hullBuilding.push_back(tri);
533  myVec3 b3(xCoords[*thirdNode], yCoords[*thirdNode], zCoords[*thirdNode]);
534  aggNodes.erase(thirdNode);
535  //Find the fourth node (most distant from triangle) to form tetrahedron
536  maxDist = 0;
537  int fourthVertex = -1;
538  for(Iter node = aggNodes.begin(); node != aggNodes.end(); node++) {
539  myVec3 thisNode(xCoords[*node], yCoords[*node], zCoords[*node]);
540  double nodeDist = pointDistFromTri(thisNode, b1, b2, b3);
541  if(nodeDist > maxDist) {
542  maxDist = nodeDist;
543  fourthVertex = *node;
544  }
545  }
546  aggNodes.remove(fourthVertex);
547  myVec3 b4(xCoords[fourthVertex], yCoords[fourthVertex], zCoords[fourthVertex]);
548  //Add three new triangles to hullBuilding to form the first tetrahedron
549  //use tri to hold the triangle info temporarily before being added to list
550  tri = hullBuilding.front();
551  tri.v1 = fourthVertex;
552  hullBuilding.push_back(tri);
553  tri = hullBuilding.front();
554  tri.v2 = fourthVertex;
555  hullBuilding.push_back(tri);
556  tri = hullBuilding.front();
557  tri.v3 = fourthVertex;
558  hullBuilding.push_back(tri);
559  //now orient all four triangles so that the vertices are oriented clockwise (so getNorm_ points outward for each)
560  myVec3 barycenter((b1.x + b2.x + b3.x + b4.x) / 4.0, (b1.y + b2.y + b3.y + b4.y) / 4.0, (b1.z + b2.z + b3.z + b4.z) / 4.0);
561  for(std::list<myTriangle>::iterator tetTri = hullBuilding.begin(); tetTri != hullBuilding.end(); tetTri++) {
562  myVec3 triVert1(xCoords[tetTri->v1], yCoords[tetTri->v1], zCoords[tetTri->v1]);
563  myVec3 triVert2(xCoords[tetTri->v2], yCoords[tetTri->v2], zCoords[tetTri->v2]);
564  myVec3 triVert3(xCoords[tetTri->v3], yCoords[tetTri->v3], zCoords[tetTri->v3]);
565  myVec3 trinorm = getNorm(triVert1, triVert2, triVert3);
566  myVec3 ptInPlane = tetTri == hullBuilding.begin() ? b1 : b4; //first triangle definitely has b1 in it, other three definitely have b4
567  if(isInFront(barycenter, ptInPlane, trinorm)) {
568  //don't want the faces of the tetrahedron to face inwards (towards barycenter) so reverse orientation
569  //by swapping two vertices
570  int temp = tetTri->v1;
571  tetTri->v1 = tetTri->v2;
572  tetTri->v2 = temp;
573  }
574  }
575  //now, have starting polyhedron in hullBuilding (all faces are facing outwards according to getNorm_) and remaining nodes to process are in aggNodes
576  //recursively, for each triangle, make a list of the points that are 'in front' of the triangle. Find the point with the maximum distance from the triangle.
577  //Add three new triangles, each sharing one edge with the original triangle but now with the most distant point as a vertex. Remove the most distant point from
578  //the list of all points that need to be processed. Also from that list remove all points that are in front of the original triangle but not in front of all three
579  //new triangles, since they are now enclosed in the hull.
580  //Construct point lists for each face of the tetrahedron individually.
581  myVec3 trinorms[4]; //normals to the four tetrahedron faces, now oriented outwards
582  int index = 0;
583  for(std::list<myTriangle>::iterator tetTri = hullBuilding.begin(); tetTri != hullBuilding.end(); tetTri++) {
584  myVec3 triVert1(xCoords[tetTri->v1], yCoords[tetTri->v1], zCoords[tetTri->v1]);
585  myVec3 triVert2(xCoords[tetTri->v2], yCoords[tetTri->v2], zCoords[tetTri->v2]);
586  myVec3 triVert3(xCoords[tetTri->v3], yCoords[tetTri->v3], zCoords[tetTri->v3]);
587  trinorms[index] = getNorm(triVert1, triVert2, triVert3);
588  index++;
589  }
590  std::list<int> startPoints1;
591  std::list<int> startPoints2;
592  std::list<int> startPoints3;
593  std::list<int> startPoints4;
594  //scope this so that 'point' is not in function scope
595  {
596  Iter point = aggNodes.begin();
597  while(!aggNodes.empty()) //this removes points one at a time as they are put in startPointsN or are already done
598  {
599  myVec3 pointVec(xCoords[*point], yCoords[*point], zCoords[*point]);
600  //Note: Because of the way the tetrahedron faces are constructed above,
601  //face 1 definitely contains b1 and faces 2-4 definitely contain b4.
602  if(isInFront(pointVec, b1, trinorms[0])) {
603  startPoints1.push_back(*point);
604  point = aggNodes.erase(point);
605  } else if(isInFront(pointVec, b4, trinorms[1])) {
606  startPoints2.push_back(*point);
607  point = aggNodes.erase(point);
608  } else if(isInFront(pointVec, b4, trinorms[2])) {
609  startPoints3.push_back(*point);
610  point = aggNodes.erase(point);
611  } else if(isInFront(pointVec, b4, trinorms[3])) {
612  startPoints4.push_back(*point);
613  point = aggNodes.erase(point);
614  } else {
615  point = aggNodes.erase(point); //points here are already inside tetrahedron.
616  }
617  }
618  //Call processTriangle for each triangle in the initial tetrahedron, one at a time.
619  }
620  typedef std::list<myTriangle>::iterator TriIter;
621  TriIter firstTri = hullBuilding.begin();
622  myTriangle start1 = *firstTri;
623  firstTri++;
624  myTriangle start2 = *firstTri;
625  firstTri++;
626  myTriangle start3 = *firstTri;
627  firstTri++;
628  myTriangle start4 = *firstTri;
629  //kick off depth-first recursive filling of hullBuilding list with all triangles in the convex hull
630  if(!startPoints1.empty())
631  processTriangle(hullBuilding, start1, startPoints1, barycenter, xCoords, yCoords, zCoords);
632  if(!startPoints2.empty())
633  processTriangle(hullBuilding, start2, startPoints2, barycenter, xCoords, yCoords, zCoords);
634  if(!startPoints3.empty())
635  processTriangle(hullBuilding, start3, startPoints3, barycenter, xCoords, yCoords, zCoords);
636  if(!startPoints4.empty())
637  processTriangle(hullBuilding, start4, startPoints4, barycenter, xCoords, yCoords, zCoords);
638  //hullBuilding now has all triangles that make up this hull.
639  //Dump hullBuilding info into the list of all triangles for the scene.
640  vertices.reserve(vertices.size() + 3 * hullBuilding.size());
641  for(TriIter hullTri = hullBuilding.begin(); hullTri != hullBuilding.end(); hullTri++) {
642  vertices.push_back(hullTri->v1);
643  vertices.push_back(hullTri->v2);
644  vertices.push_back(hullTri->v3);
645  geomSizes.push_back(3);
646  }
647  }
648  }
649 
650 #ifdef HAVE_MUELU_CGAL
651  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
652  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doCGALConvexHulls3D(std::vector<int>& vertices, std::vector<int>& geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector<bool>& isRoot, const ArrayRCP<LO>& vertex2AggId, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& xCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& yCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& zCoords) {
653 
654  typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
655  typedef K::Point_3 Point_3;
656  typedef CGAL::Polyhedron_3<K> Polyhedron_3;
657  typedef std::vector<int>::iterator Iter;
658  for(int agg = 0; agg < numLocalAggs; agg++) {
659  std::vector<int> aggNodes;
660  std::vector<Point_3> aggPoints;
661  for(int i = 0; i < numFineNodes; i++) {
662  if(vertex2AggId[i] == agg) {
663  Point_3 p(xCoords[i], yCoords[i], zCoords[i]);
664  aggPoints.push_back(p);
665  aggNodes.push_back(i);
666  }
667  }
668  //First, check anomalous cases
670  "CoarseningVisualization::doCGALConvexHulls3D: aggregate contains zero nodes!");
671  if(aggNodes.size() == 1) {
672  vertices.push_back(aggNodes.front());
673  geomSizes.push_back(1);
674  continue;
675  } else if(aggNodes.size() == 2) {
676  vertices.push_back(aggNodes.front());
677  vertices.push_back(aggNodes.back());
678  geomSizes.push_back(2);
679  continue;
680  }
681  //check for collinearity
682  bool areCollinear = true;
683  {
684  Iter it = aggNodes.begin();
685  myVec3 firstVec(xCoords[*it], yCoords[*it], zCoords[*it]);
686  myVec3 comp;
687  {
688  it++;
689  myVec3 secondVec(xCoords[*it], yCoords[*it], zCoords[*it]); //cross this with other vectors to compare
690  comp = vecSubtract(secondVec, firstVec);
691  it++;
692  }
693  for(; it != aggNodes.end(); it++) {
694  myVec3 thisVec(xCoords[*it], yCoords[*it], zCoords[*it]);
695  myVec3 cross = crossProduct(vecSubtract(thisVec, firstVec), comp);
696  if(mymagnitude(cross) > 1e-8) {
697  areCollinear = false;
698  break;
699  }
700  }
701  }
702  if(areCollinear)
703  {
704  //find the endpoints of segment describing all the points
705  //compare x, if tie compare y, if tie compare z
706  Iter min = aggNodes.begin();
707  Iter max = aggNodes.begin();
708  Iter it = ++aggNodes.begin();
709  for(; it != aggNodes.end(); it++) {
710  if(xCoords[*it] < xCoords[*min]) min = it;
711  else if(xCoords[*it] == xCoords[*min]) {
712  if(yCoords[*it] < yCoords[*min]) min = it;
713  else if(yCoords[*it] == yCoords[*min]) {
714  if(zCoords[*it] < zCoords[*min]) min = it;
715  }
716  }
717  if(xCoords[*it] > xCoords[*max]) max = it;
718  else if(xCoords[*it] == xCoords[*max]) {
719  if(yCoords[*it] > yCoords[*max]) max = it;
720  else if(yCoords[*it] == yCoords[*max]) {
721  if(zCoords[*it] > zCoords[*max])
722  max = it;
723  }
724  }
725  }
726  vertices.push_back(*min);
727  vertices.push_back(*max);
728  geomSizes.push_back(2);
729  continue;
730  }
731  // do not handle coplanar or general case here. Just let's use CGAL
732  {
733  Polyhedron_3 result;
734  CGAL::convex_hull_3( aggPoints.begin(), aggPoints.end(), result);
735 
736  // loop over all facets
737  Polyhedron_3::Facet_iterator fi;
738  for (fi = result.facets_begin(); fi != result.facets_end(); fi++) {
739  int cntVertInAgg = 0;
740  Polyhedron_3::Halfedge_around_facet_const_circulator hit = fi->facet_begin();
741  do {
742  const Point_3 & pp = hit->vertex()->point();
743  // loop over all aggregate nodes and find corresponding node id
744  for(size_t l = 0; l < aggNodes.size(); l++)
745  {
746  if(fabs(pp.x() - xCoords[aggNodes[l]]) < 1e-12 &&
747  fabs(pp.y() - yCoords[aggNodes[l]]) < 1e-12 &&
748  fabs(pp.z() - zCoords[aggNodes[l]]) < 1e-12)
749  {
750  vertices.push_back(aggNodes[l]);
751  cntVertInAgg++;
752  break;
753  }
754  }
755  } while (++hit != fi->facet_begin());
756  geomSizes.push_back(cntVertInAgg);
757  }
758  }
759  }
760 
761  }
762 #endif
763 
764  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
768 
769  std::vector<std::pair<int, int> > vert1; //vertices (node indices)
770 
772  for(LocalOrdinal locRow = 0; locRow < LocalOrdinal(G->GetNodeNumVertices()); locRow++) {
773  neighbors = G->getNeighborVertices(locRow);
774  //Add those local indices (columns) to the list of connections (which are pairs of the form (localM, localN))
775  for(int gEdge = 0; gEdge < int(neighbors.size()); ++gEdge) {
776  vert1.push_back(std::pair<int, int>(locRow, neighbors[gEdge]));
777  }
778  }
779  for(size_t i = 0; i < vert1.size(); i ++) {
780  if(vert1[i].first > vert1[i].second) {
781  int temp = vert1[i].first;
782  vert1[i].first = vert1[i].second;
783  vert1[i].second = temp;
784  }
785  }
786  std::sort(vert1.begin(), vert1.end());
787  std::vector<std::pair<int, int> >::iterator newEnd = unique(vert1.begin(), vert1.end()); //remove duplicate edges
788  vert1.erase(newEnd, vert1.end());
789  //std::vector<int> points1;
790  vertices.reserve(2 * vert1.size());
791  geomSizes.reserve(vert1.size());
792  for(size_t i = 0; i < vert1.size(); i++) {
793  vertices.push_back(vert1[i].first);
794  vertices.push_back(vert1[i].second);
795  geomSizes.push_back(2);
796  }
797  }
798 
799  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
801  {
802  return myVec3(v1.y * v2.z - v1.z * v2.y, v1.z * v2.x - v1.x * v2.z, v1.x * v2.y - v1.y * v2.x);
803  }
804 
805  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
807  {
808  return v1.x * v2.x + v1.y * v2.y;
809  }
810 
811  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
813  {
814  return v1.x * v2.x + v1.y * v2.y + v1.z * v2.z;
815  }
816 
817  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
819  {
820  myVec3 rel(point.x - inPlane.x, point.y - inPlane.y, point.z - inPlane.z); //position of the point relative to the plane
821  return dotProduct(rel, n) > 1e-12 ? true : false;
822  }
823 
824  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
826  {
827  return sqrt(dotProduct(vec, vec));
828  }
829 
830  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
832  {
833  return sqrt(dotProduct(vec, vec));
834  }
835 
836  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
838  {
839  return mymagnitude(vecSubtract(p1, p2));
840  }
841 
842 
843  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
845  {
846  return mymagnitude(vecSubtract(p1, p2));
847  }
848 
849  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
851  {
852  return myVec2(v1.x - v2.x, v1.y - v2.y);
853  }
854 
855  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
857  {
858  return myVec3(v1.x - v2.x, v1.y - v2.y, v1.z - v2.z);
859  }
860 
861  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
862  myVec2 VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNorm(myVec2 v) //"normal" to a 2D vector - just rotate 90 degrees to left
863  {
864  return myVec2(v.y, -v.x);
865  }
866 
867  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
868  myVec3 VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNorm(myVec3 v1, myVec3 v2, myVec3 v3) //normal to face of triangle (will be outward rel. to polyhedron) (v1, v2, v3 are in CCW order when normal is toward viewpoint)
869  {
870  return crossProduct(vecSubtract(v2, v1), vecSubtract(v3, v1));
871  }
872 
873  //get minimum distance from 'point' to plane containing v1, v2, v3 (or the triangle with v1, v2, v3 as vertices)
874  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
876  {
877  using namespace std;
878  myVec3 norm = getNorm(v1, v2, v3);
879  //must normalize the normal vector
880  double normScl = mymagnitude(norm);
881  double rv = 0.0;
882  if (normScl > 1e-8) {
883  norm.x /= normScl;
884  norm.y /= normScl;
885  norm.z /= normScl;
886  rv = fabs(dotProduct(norm, vecSubtract(point, v1)));
887  } else {
888  // triangle is degenerated
889  myVec3 test1 = vecSubtract(v3, v1);
890  myVec3 test2 = vecSubtract(v2, v1);
891  bool useTest1 = mymagnitude(test1) > 0.0 ? true : false;
892  bool useTest2 = mymagnitude(test2) > 0.0 ? true : false;
893  if(useTest1 == true) {
894  double normScl1 = mymagnitude(test1);
895  test1.x /= normScl1;
896  test1.y /= normScl1;
897  test1.z /= normScl1;
898  rv = fabs(dotProduct(test1, vecSubtract(point, v1)));
899  } else if (useTest2 == true) {
900  double normScl2 = mymagnitude(test2);
901  test2.x /= normScl2;
902  test2.y /= normScl2;
903  test2.z /= normScl2;
904  rv = fabs(dotProduct(test2, vecSubtract(point, v1)));
905  } else {
907  "VisualizationHelpers::pointDistFromTri: Could not determine the distance of a point to a triangle.");
908  }
909 
910  }
911  return rv;
912  }
913 
914  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
915  std::vector<myTriangle> VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::processTriangle(std::list<myTriangle>& tris, myTriangle tri, std::list<int>& pointsInFront, myVec3& barycenter, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& xCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& yCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>& zCoords) {
916  //*tri is in the tris list, and is the triangle to process here. tris is a complete list of all triangles in the hull so far. pointsInFront is only a list of the nodes in front of tri. Need coords also.
917  //precondition: each triangle is already oriented so that getNorm_(v1, v2, v3) points outward (away from interior of hull)
918  //First find the point furthest from triangle.
919  using namespace std;
920  typedef std::list<int>::iterator Iter;
921  typedef std::list<myTriangle>::iterator TriIter;
922  typedef list<pair<int, int> >::iterator EdgeIter;
923  double maxDist = 0;
924  //Need vector representations of triangle's vertices
925  myVec3 v1(xCoords[tri.v1], yCoords[tri.v1], zCoords[tri.v1]);
926  myVec3 v2(xCoords[tri.v2], yCoords[tri.v2], zCoords[tri.v2]);
927  myVec3 v3(xCoords[tri.v3], yCoords[tri.v3], zCoords[tri.v3]);
928  myVec3 farPointVec; //useful to have both the point's coordinates and it's position in the list
929  Iter farPoint = pointsInFront.begin();
930  for(Iter point = pointsInFront.begin(); point != pointsInFront.end(); point++)
931  {
932  myVec3 pointVec(xCoords[*point], yCoords[*point], zCoords[*point]);
933  double dist = pointDistFromTri(pointVec, v1, v2, v3);
934  if(dist > maxDist)
935  {
936  dist = maxDist;
937  farPointVec = pointVec;
938  farPoint = point;
939  }
940  }
941  //Find all the triangles that the point is in front of (can be more than 1)
942  //At the same time, remove them from tris, as every one will be replaced later
943  vector<myTriangle> visible; //use a list of iterators so that the underlying object is still in tris
944  for(TriIter it = tris.begin(); it != tris.end();)
945  {
946  myVec3 vec1(xCoords[it->v1], yCoords[it->v1], zCoords[it->v1]);
947  myVec3 vec2(xCoords[it->v2], yCoords[it->v2], zCoords[it->v2]);
948  myVec3 vec3(xCoords[it->v3], yCoords[it->v3], zCoords[it->v3]);
949  myVec3 norm = getNorm(vec1, vec2, vec3);
950  if(isInFront(farPointVec, vec1, norm))
951  {
952  visible.push_back(*it);
953  it = tris.erase(it);
954  }
955  else
956  it++;
957  }
958  //Figure out what triangles need to be destroyed/created
959  //First create a list of edges (as std::pair<int, int>, where the two ints are the node endpoints)
960  list<pair<int, int> > horizon;
961  //For each triangle, add edges to the list iff the edge only appears once in the set of all
962  //Have members of horizon have the lower node # first, and the higher one second
963  for(vector<myTriangle>::iterator it = visible.begin(); it != visible.end(); it++)
964  {
965  pair<int, int> e1(it->v1, it->v2);
966  pair<int, int> e2(it->v2, it->v3);
967  pair<int, int> e3(it->v1, it->v3);
968  //"sort" the pair values
969  if(e1.first > e1.second)
970  {
971  int temp = e1.first;
972  e1.first = e1.second;
973  e1.second = temp;
974  }
975  if(e2.first > e2.second)
976  {
977  int temp = e2.first;
978  e2.first = e2.second;
979  e2.second = temp;
980  }
981  if(e3.first > e3.second)
982  {
983  int temp = e3.first;
984  e3.first = e3.second;
985  e3.second = temp;
986  }
987  horizon.push_back(e1);
988  horizon.push_back(e2);
989  horizon.push_back(e3);
990  }
991  //sort based on lower node first, then higher node (lexicographical ordering built in to pair)
992  horizon.sort();
993  //Remove all edges from horizon, except those that appear exactly once
994  {
995  EdgeIter it = horizon.begin();
996  while(it != horizon.end())
997  {
998  int occur = count(horizon.begin(), horizon.end(), *it);
999  if(occur > 1)
1000  {
1001  pair<int, int> removeVal = *it;
1002  while(removeVal == *it && !(it == horizon.end()))
1003  it = horizon.erase(it);
1004  }
1005  else
1006  it++;
1007  }
1008  }
1009  //Now make a list of new triangles being created, each of which take 2 vertices from an edge and one from farPoint
1010  list<myTriangle> newTris;
1011  for(EdgeIter it = horizon.begin(); it != horizon.end(); it++)
1012  {
1013  myTriangle t(it->first, it->second, *farPoint);
1014  newTris.push_back(t);
1015  }
1016  //Ensure every new triangle is oriented outwards, using the barycenter of the initial tetrahedron
1017  vector<myTriangle> trisToProcess;
1018  vector<list<int> > newFrontPoints;
1019  for(TriIter it = newTris.begin(); it != newTris.end(); it++)
1020  {
1021  myVec3 t1(xCoords[it->v1], yCoords[it->v1], zCoords[it->v1]);
1022  myVec3 t2(xCoords[it->v2], yCoords[it->v2], zCoords[it->v2]);
1023  myVec3 t3(xCoords[it->v3], yCoords[it->v3], zCoords[it->v3]);
1024  if(isInFront(barycenter, t1, getNorm(t1, t2, t3)))
1025  {
1026  //need to swap two vertices to flip orientation of triangle
1027  int temp = it->v1;
1028  myVec3 tempVec = t1;
1029  it->v1 = it->v2;
1030  t1 = t2;
1031  it->v2 = temp;
1032  t2 = tempVec;
1033  }
1034  myVec3 outwardNorm = getNorm(t1, t2, t3); //now definitely points outwards
1035  //Add the triangle to tris
1036  tris.push_back(*it);
1037  trisToProcess.push_back(tris.back());
1038  //Make a list of the points that are in front of nextToProcess, to be passed in for processing
1039  list<int> newInFront;
1040  for(Iter point = pointsInFront.begin(); point != pointsInFront.end();)
1041  {
1042  myVec3 pointVec(xCoords[*point], yCoords[*point], zCoords[*point]);
1043  if(isInFront(pointVec, t1, outwardNorm))
1044  {
1045  newInFront.push_back(*point);
1046  point = pointsInFront.erase(point);
1047  }
1048  else
1049  point++;
1050  }
1051  newFrontPoints.push_back(newInFront);
1052  }
1053  vector<myTriangle> allRemoved; //list of all invalid iterators that were erased by calls to processmyTriangle below
1054  for(int i = 0; i < int(trisToProcess.size()); i++)
1055  {
1056  if(!newFrontPoints[i].empty())
1057  {
1058  //Comparing the 'triangle to process' to the one for this call prevents infinite recursion/stack overflow.
1059  //TODO: Why was it doing that? Rounding error? Make more robust fix. But this does work for the time being.
1060  if(find(allRemoved.begin(), allRemoved.end(), trisToProcess[i]) == allRemoved.end() && !(trisToProcess[i] == tri))
1061  {
1062  vector<myTriangle> removedList = processTriangle(tris, trisToProcess[i], newFrontPoints[i], barycenter, xCoords, yCoords, zCoords);
1063  for(int j = 0; j < int(removedList.size()); j++)
1064  allRemoved.push_back(removedList[j]);
1065  }
1066  }
1067  }
1068  return visible;
1069  }
1070 
1071  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1072  std::vector<int> VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::giftWrap(std::vector<myVec2>& points, std::vector<int>& nodes, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType> & xCoords, const Teuchos::ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType> & yCoords) {
1074  "CoarseningVisualization::giftWrap: Gift wrap algorithm input has to have at least 3 points!");
1075 
1076 #if 1 // TAW's version to determine "minimal" node
1077  // determine minimal x and y coordinates
1078  double min_x =points[0].x;
1079  double min_y =points[0].y;
1080  for(std::vector<int>::iterator it = nodes.begin(); it != nodes.end(); it++) {
1081  int i = it - nodes.begin();
1082  if(points[i].x < min_x) min_x = points[i].x;
1083  if(points[i].y < min_y) min_y = points[i].y;
1084  }
1085  // create dummy min coordinates
1086  min_x -= 1.0;
1087  min_y -= 1.0;
1088  myVec2 dummy_min(min_x, min_y);
1089 
1090  // loop over all nodes and determine nodes with minimal distance to (min_x, min_y)
1091  std::vector<int> hull;
1092  myVec2 min = points[0];
1093  double mindist = distance(min,dummy_min);
1094  std::vector<int>::iterator minNode = nodes.begin();
1095  for(std::vector<int>::iterator it = nodes.begin(); it != nodes.end(); it++) {
1096  int i = it - nodes.begin();
1097  if(distance(points[i],dummy_min) < mindist) {
1098  mindist = distance(points[i],dummy_min);
1099  min = points[i];
1100  minNode = it;
1101  }
1102  }
1103  hull.push_back(*minNode);
1104 #else // Brian's code
1105  std::vector<int> hull;
1106  std::vector<int>::iterator minNode = nodes.begin();
1107  myVec2 min = points[0];
1108  for(std::vector<int>::iterator it = nodes.begin(); it != nodes.end(); it++)
1109  {
1110  int i = it - nodes.begin();
1111  if(points[i].x < min.x || (fabs(points[i].x - min.x) < 1e-10 && points[i].y < min.y))
1112  {
1113  min = points[i];
1114  minNode = it;
1115  }
1116  }
1117  hull.push_back(*minNode);
1118 #endif
1119 
1120  bool includeMin = false;
1121  //int debug_it = 0;
1122  while(1)
1123  {
1124  std::vector<int>::iterator leftMost = nodes.begin();
1125  if(!includeMin && leftMost == minNode)
1126  {
1127  leftMost++;
1128  }
1129  std::vector<int>::iterator it = leftMost;
1130  it++;
1131  for(; it != nodes.end(); it++)
1132  {
1133  if(it == minNode && !includeMin) //don't compare to min on very first sweep
1134  continue;
1135  if(*it == hull.back())
1136  continue;
1137  //see if it is in front of line containing nodes thisHull.back() and leftMost
1138  //first get the left normal of leftMost - thisHull.back() (<dy, -dx>)
1139  myVec2 leftMostVec = points[leftMost - nodes.begin()];
1140  myVec2 lastVec(xCoords[hull.back()], yCoords[hull.back()]);
1141  myVec2 testNorm = getNorm(vecSubtract(leftMostVec, lastVec));
1142  //now dot testNorm with *it - leftMost. If dot is positive, leftMost becomes it. If dot is zero, take one further from thisHull.back().
1143  myVec2 itVec(xCoords[*it], yCoords[*it]);
1144  double dotProd = dotProduct(testNorm, vecSubtract(itVec, lastVec));
1145  if(-1e-8 < dotProd && dotProd < 1e-8)
1146  {
1147  //thisHull.back(), it and leftMost are collinear.
1148  //Just sum the differences in x and differences in y for each and compare to get further one, don't need distance formula
1149  myVec2 itDist = vecSubtract(itVec, lastVec);
1150  myVec2 leftMostDist = vecSubtract(leftMostVec, lastVec);
1151  if(fabs(itDist.x) + fabs(itDist.y) > fabs(leftMostDist.x) + fabs(leftMostDist.y)) {
1152  leftMost = it;
1153  }
1154  }
1155  else if(dotProd > 0) {
1156  leftMost = it;
1157 
1158  }
1159  }
1160  //if leftMost is min, then the loop is complete.
1161  if(*leftMost == *minNode)
1162  break;
1163  hull.push_back(*leftMost);
1164  includeMin = true; //have found second point (the one after min) so now include min in the searches
1165  //debug_it ++;
1166  //if(debug_it > 100) exit(0); //break;
1167  }
1168  return hull;
1169  }
1170 
1171  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1172  std::vector<int> VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::makeUnique(std::vector<int>& vertices) const
1173  {
1174  using namespace std;
1175  vector<int> uniqueNodes = vertices;
1176  sort(uniqueNodes.begin(), uniqueNodes.end());
1177  vector<int>::iterator newUniqueFineEnd = unique(uniqueNodes.begin(), uniqueNodes.end());
1178  uniqueNodes.erase(newUniqueFineEnd, uniqueNodes.end());
1179  //uniqueNodes is now a sorted list of the nodes whose info actually goes in file
1180  //Now replace values in vertices with locations of the old values in uniqueFine
1181  for(int i = 0; i < int(vertices.size()); i++)
1182  {
1183  int lo = 0;
1184  int hi = uniqueNodes.size() - 1;
1185  int mid = 0;
1186  int search = vertices[i];
1187  while(lo <= hi)
1188  {
1189  mid = lo + (hi - lo) / 2;
1190  if(uniqueNodes[mid] == search)
1191  break;
1192  else if(uniqueNodes[mid] > search)
1193  hi = mid - 1;
1194  else
1195  lo = mid + 1;
1196  }
1197  if(uniqueNodes[mid] != search)
1198  throw runtime_error("Issue in makeUnique_() - a point wasn't found in list.");
1199  vertices[i] = mid;
1200  }
1201  return uniqueNodes;
1202  }
1203 
1204  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1205  std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceAll(std::string result, const std::string& replaceWhat, const std::string& replaceWithWhat) const {
1206  while(1) {
1207  const int pos = result.find(replaceWhat);
1208  if (pos == -1)
1209  break;
1210  result.replace(pos, replaceWhat.size(), replaceWithWhat);
1211  }
1212  return result;
1213  }
1214 
1215  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1216  std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getFileName(int numProcs, int myRank, int level, const Teuchos::ParameterList & pL) const {
1217  std::string filenameToWrite = getBaseFileName(numProcs, level, pL);
1218  filenameToWrite = this->replaceAll(filenameToWrite, "%PROCID", toString(myRank));
1219  return filenameToWrite;
1220  }
1221 
1222  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1224  std::string filenameToWrite = pL.get<std::string>("visualization: output filename");
1225  int timeStep = pL.get<int>("visualization: output file: time step");
1226  int iter = pL.get<int>("visualization: output file: iter");
1227 
1228  if(filenameToWrite.rfind(".vtu") == std::string::npos)
1229  filenameToWrite.append(".vtu");
1230  if(numProcs > 1 && filenameToWrite.rfind("%PROCID") == std::string::npos) //filename can't be identical between processsors in parallel problem
1231  filenameToWrite.insert(filenameToWrite.rfind(".vtu"), "-proc%PROCID");
1232 
1233  filenameToWrite = this->replaceAll(filenameToWrite, "%LEVELID", toString(level));
1234  filenameToWrite = this->replaceAll(filenameToWrite, "%TIMESTEP", toString(timeStep));
1235  filenameToWrite = this->replaceAll(filenameToWrite, "%ITER", toString(iter));
1236  return filenameToWrite;
1237  }
1238 
1239  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1240  std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getPVTUFileName(int numProcs, int /* myRank */, int level, const Teuchos::ParameterList & pL) const {
1241  std::string filenameToWrite = getBaseFileName(numProcs, level, pL);
1242  std::string masterStem = filenameToWrite.substr(0, filenameToWrite.rfind(".vtu"));
1243  masterStem = this->replaceAll(masterStem, "%PROCID", "");
1244  std::string pvtuFilename = masterStem + "-master.pvtu";
1245  return pvtuFilename;
1246  }
1247 
1248  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1249  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKOpening(std::ofstream & fout, std::vector<int> & uniqueFine, std::vector<int> & geomSizesFine) const {
1250  std::string styleName = "PointCloud"; // TODO adapt this
1251 
1252  std::string indent = " ";
1253  fout << "<!--" << styleName << " Aggregates Visualization-->" << std::endl;
1254  fout << "<VTKFile type=\"UnstructuredGrid\" byte_order=\"LittleEndian\">" << std::endl;
1255  fout << " <UnstructuredGrid>" << std::endl;
1256  fout << " <Piece NumberOfPoints=\"" << uniqueFine.size() << "\" NumberOfCells=\"" << geomSizesFine.size() << "\">" << std::endl;
1257  fout << " <PointData Scalars=\"Node Aggregate Processor\">" << std::endl;
1258  }
1259 
1260  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1261  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKNodes(std::ofstream & fout, std::vector<int> & uniqueFine, Teuchos::RCP<const Map> & nodeMap) const {
1262  std::string indent = " ";
1263  fout << " <DataArray type=\"Int32\" Name=\"Node\" format=\"ascii\">" << std::endl;
1264  indent = " ";
1265  bool localIsGlobal = GlobalOrdinal(nodeMap->getGlobalNumElements()) == GlobalOrdinal(nodeMap->getNodeNumElements());
1266  for(size_t i = 0; i < uniqueFine.size(); i++)
1267  {
1268  if(localIsGlobal)
1269  {
1270  fout << uniqueFine[i] << " "; //if all nodes are on this processor, do not map from local to global
1271  }
1272  else
1273  fout << nodeMap->getGlobalElement(uniqueFine[i]) << " ";
1274  if(i % 10 == 9)
1275  fout << std::endl << indent;
1276  }
1277  fout << std::endl;
1278  fout << " </DataArray>" << std::endl;
1279  }
1280 
1281  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1282  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKData(std::ofstream & fout, std::vector<int> & uniqueFine, LocalOrdinal myAggOffset, ArrayRCP<LocalOrdinal> & vertex2AggId, int myRank) const {
1283  std::string indent = " ";
1284  fout << " <DataArray type=\"Int32\" Name=\"Aggregate\" format=\"ascii\">" << std::endl;
1285  fout << indent;
1286  for(int i = 0; i < int(uniqueFine.size()); i++)
1287  {
1288  fout << myAggOffset + vertex2AggId[uniqueFine[i]] << " ";
1289  if(i % 10 == 9)
1290  fout << std::endl << indent;
1291  }
1292  fout << std::endl;
1293  fout << " </DataArray>" << std::endl;
1294  fout << " <DataArray type=\"Int32\" Name=\"Processor\" format=\"ascii\">" << std::endl;
1295  fout << indent;
1296  for(int i = 0; i < int(uniqueFine.size()); i++)
1297  {
1298  fout << myRank << " ";
1299  if(i % 20 == 19)
1300  fout << std::endl << indent;
1301  }
1302  fout << std::endl;
1303  fout << " </DataArray>" << std::endl;
1304  fout << " </PointData>" << std::endl;
1305  }
1306 
1307  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1309  std::string indent = " ";
1310  fout << " <Points>" << std::endl;
1311  fout << " <DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">" << std::endl;
1312  fout << indent;
1313  for(int i = 0; i < int(uniqueFine.size()); i++)
1314  {
1315  fout << fx[uniqueFine[i]] << " " << fy[uniqueFine[i]] << " ";
1316  if(dim == 2)
1317  fout << "0 ";
1318  else
1319  fout << fz[uniqueFine[i]] << " ";
1320  if(i % 3 == 2)
1321  fout << std::endl << indent;
1322  }
1323  fout << std::endl;
1324  fout << " </DataArray>" << std::endl;
1325  fout << " </Points>" << std::endl;
1326  }
1327 
1328  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1329  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKCells(std::ofstream & fout, std::vector<int> & /* uniqueFine */, std::vector<LocalOrdinal> & vertices, std::vector<LocalOrdinal> & geomSize) const {
1330  std::string indent = " ";
1331  fout << " <Cells>" << std::endl;
1332  fout << " <DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">" << std::endl;
1333  fout << indent;
1334  for(int i = 0; i < int(vertices.size()); i++)
1335  {
1336  fout << vertices[i] << " ";
1337  if(i % 10 == 9)
1338  fout << std::endl << indent;
1339  }
1340  fout << std::endl;
1341  fout << " </DataArray>" << std::endl;
1342  fout << " <DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">" << std::endl;
1343  fout << indent;
1344  int accum = 0;
1345  for(int i = 0; i < int(geomSize.size()); i++)
1346  {
1347  accum += geomSize[i];
1348  fout << accum << " ";
1349  if(i % 10 == 9)
1350  fout << std::endl << indent;
1351  }
1352  fout << std::endl;
1353  fout << " </DataArray>" << std::endl;
1354  fout << " <DataArray type=\"Int32\" Name=\"types\" format=\"ascii\">" << std::endl;
1355  fout << indent;
1356  for(int i = 0; i < int(geomSize.size()); i++)
1357  {
1358  switch(geomSize[i])
1359  {
1360  case 1:
1361  fout << "1 "; //Point
1362  break;
1363  case 2:
1364  fout << "3 "; //Line
1365  break;
1366  case 3:
1367  fout << "5 "; //Triangle
1368  break;
1369  default:
1370  fout << "7 "; //Polygon
1371  }
1372  if(i % 30 == 29)
1373  fout << std::endl << indent;
1374  }
1375  fout << std::endl;
1376  fout << " </DataArray>" << std::endl;
1377  fout << " </Cells>" << std::endl;
1378  }
1379 
1380  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1382  fout << " </Piece>" << std::endl;
1383  fout << " </UnstructuredGrid>" << std::endl;
1384  fout << "</VTKFile>" << std::endl;
1385  }
1386 
1387 
1388  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1389  void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writePVTU(std::ofstream& pvtu, std::string baseFname, int numProcs, bool bFineEdges, bool /* bCoarseEdges */) const {
1390  //If using vtk, filenameToWrite now contains final, correct ***.vtu filename (for the current proc)
1391  //So the root proc will need to use its own filenameToWrite to make a list of the filenames of all other procs to put in
1392  //pvtu file.
1393  pvtu << "<VTKFile type=\"PUnstructuredGrid\" byte_order=\"LittleEndian\">" << std::endl;
1394  pvtu << " <PUnstructuredGrid GhostLevel=\"0\">" << std::endl;
1395  pvtu << " <PPointData Scalars=\"Node Aggregate Processor\">" << std::endl;
1396  pvtu << " <PDataArray type=\"Int32\" Name=\"Node\"/>" << std::endl;
1397  pvtu << " <PDataArray type=\"Int32\" Name=\"Aggregate\"/>" << std::endl;
1398  pvtu << " <PDataArray type=\"Int32\" Name=\"Processor\"/>" << std::endl;
1399  pvtu << " </PPointData>" << std::endl;
1400  pvtu << " <PPoints>" << std::endl;
1401  pvtu << " <PDataArray type=\"Float64\" NumberOfComponents=\"3\"/>" << std::endl;
1402  pvtu << " </PPoints>" << std::endl;
1403  for(int i = 0; i < numProcs; i++) {
1404  //specify the piece for each proc (the replaceAll expression matches what the filenames will be on other procs)
1405  pvtu << " <Piece Source=\"" << replaceAll(baseFname, "%PROCID", toString(i)) << "\"/>" << std::endl;
1406  }
1407  //reference files with graph pieces, if applicable
1408  if(bFineEdges)
1409  {
1410  for(int i = 0; i < numProcs; i++)
1411  {
1412  std::string fn = replaceAll(baseFname, "%PROCID", toString(i));
1413  pvtu << " <Piece Source=\"" << fn.insert(fn.rfind(".vtu"), "-finegraph") << "\"/>" << std::endl;
1414  }
1415  }
1416  /*if(doCoarseGraphEdges_)
1417  {
1418  for(int i = 0; i < numProcs; i++)
1419  {
1420  std::string fn = replaceAll(baseFname, "%PROCID", toString(i));
1421  pvtu << " <Piece Source=\"" << fn.insert(fn.rfind(".vtu"), "-coarsegraph") << "\"/>" << std::endl;
1422  }
1423  }*/
1424  pvtu << " </PUnstructuredGrid>" << std::endl;
1425  pvtu << "</VTKFile>" << std::endl;
1426  pvtu.close();
1427  }
1428 
1429  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1431  try {
1432  std::ofstream color("random-colormap.xml");
1433  color << "<ColorMap name=\"MueLu-Random\" space=\"RGB\">" << std::endl;
1434  //Give -1, -2, -3 distinctive colors (so that the style functions can have constrasted geometry types)
1435  //Do red, orange, yellow to constrast with cool color spectrum for other types
1436  color << " <Point x=\"" << -1 << "\" o=\"1\" r=\"1\" g=\"0\" b=\"0\"/>" << std::endl;
1437  color << " <Point x=\"" << -2 << "\" o=\"1\" r=\"1\" g=\"0.6\" b=\"0\"/>" << std::endl;
1438  color << " <Point x=\"" << -3 << "\" o=\"1\" r=\"1\" g=\"1\" b=\"0\"/>" << std::endl;
1439  srand(time(NULL));
1440  for(int i = 0; i < 5000; i += 4) {
1441  color << " <Point x=\"" << i << "\" o=\"1\" r=\"" << (rand() % 50) / 256.0 << "\" g=\"" << (rand() % 256) / 256.0 << "\" b=\"" << (rand() % 256) / 256.0 << "\"/>" << std::endl;
1442  }
1443  color << "</ColorMap>" << std::endl;
1444  color.close();
1445  }
1446  catch(std::exception& e) {
1448  "VisualizationHelpers::buildColormap: Error while building colormap file: " << e.what());
1449  }
1450  }
1451 
1452 } // namespace MueLu
1453 
1454 #endif /* MUELU_VISUALIZATIONHELPERS_DEF_HPP_ */
static void doPointCloud(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
std::string getPVTUFileName(int numProcs, int myRank, int level, const Teuchos::ParameterList &pL) const
void writeFileVTKOpening(std::ofstream &fout, std::vector< int > &uniqueFine, std::vector< int > &geomSizesFine) const
void writeFileVTKCoordinates(std::ofstream &fout, std::vector< int > &uniqueFine, Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &fx, Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &fy, Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &fz, int dim) const
std::string replaceAll(std::string result, const std::string &replaceWhat, const std::string &replaceWithWhat) const
static void doConvexHulls2D(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector< bool > &isRoot, const ArrayRCP< LO > &vertex2AggId, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &xCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &yCoords)
static void doCGALConvexHulls3D(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector< bool > &isRoot, const ArrayRCP< LO > &vertex2AggId, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &xCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &yCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &zCoords)
T & get(const std::string &name, T def_value)
void writeFileVTKData(std::ofstream &fout, std::vector< int > &uniqueFine, LocalOrdinal myAggOffset, ArrayRCP< LocalOrdinal > &vertex2AggId, int myRank) const
static void doCGALConvexHulls2D(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector< bool > &isRoot, const ArrayRCP< LO > &vertex2AggId, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &xCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &yCoords)
virtual size_t GetNodeNumVertices() const =0
Return number of vertices owned by the calling node.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::vector< int > makeUnique(std::vector< int > &vertices) const
replaces node indices in vertices with compressed unique indices, and returns list of unique points ...
size_type size() const
std::string getFileName(int numProcs, int myRank, int level, const Teuchos::ParameterList &pL) const
LocalOrdinal LO
void writeFileVTKNodes(std::ofstream &fout, std::vector< int > &uniqueFine, Teuchos::RCP< const Map > &nodeMap) const
RCP< ParameterList > GetValidParameterList() const
static double pointDistFromTri(myVec3 point, myVec3 v1, myVec3 v2, myVec3 v3)
static double dotProduct(myVec2 v1, myVec2 v2)
void writePVTU(std::ofstream &pvtu, std::string baseFname, int numProcs, bool bFineEdges=false, bool bCoarseEdges=false) const
std::string getBaseFileName(int numProcs, int level, const Teuchos::ParameterList &pL) const
static void doConvexHulls3D(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector< bool > &isRoot, const ArrayRCP< LO > &vertex2AggId, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &xCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &yCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &zCoords)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void doGraphEdges(std::vector< int > &vertices, std::vector< int > &geomSizes, Teuchos::RCP< GraphBase > &G, Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &fx, Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &fy, Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &fz)
static void doJacks(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector< bool > &isRoot, const ArrayRCP< LO > &vertex2AggId)
static myVec3 crossProduct(myVec3 v1, myVec3 v2)
TransListIter iter
static myVec2 vecSubtract(myVec2 v1, myVec2 v2)
static double distance(myVec2 p1, myVec2 p2)
void writeFileVTKCells(std::ofstream &fout, std::vector< int > &uniqueFine, std::vector< LocalOrdinal > &vertices, std::vector< LocalOrdinal > &geomSize) const
static bool isInFront(myVec3 point, myVec3 inPlane, myVec3 n)
Exception throws to report errors in the internal logical of the program.
static std::vector< myTriangle > processTriangle(std::list< myTriangle > &tris, myTriangle tri, std::list< int > &pointsInFront, myVec3 &barycenter, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &xCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &yCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &zCoords)
void replaceAll(std::string &str, const std::string &from, const std::string &to)
virtual Teuchos::ArrayView< const LocalOrdinal > getNeighborVertices(LocalOrdinal v) const =0
Return the list of vertices adjacent to the vertex &#39;v&#39;.
static std::vector< int > giftWrap(std::vector< myVec2 > &points, std::vector< int > &nodes, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &xCoords, const Teuchos::ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::coordinateType > &yCoords)
void writeFileVTKClosing(std::ofstream &fout) const