MueLu Version of the Day
Loading...
Searching...
No Matches
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>
53#include <Xpetra_MultiVectorFactory.hpp>
55#include "MueLu_Level.hpp"
56#include "MueLu_Graph.hpp"
57#include <vector>
58#include <list>
59#include <algorithm>
60#include <string>
61
62
63namespace MueLu {
64
65 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67 RCP<ParameterList> validParamList = rcp(new ParameterList());
68
69 validParamList->set< std::string > ("visualization: output filename", "viz%LEVELID", "filename for VTK-style visualization output");
70 validParamList->set< int > ("visualization: output file: time step", 0, "time step variable for output file name");// Remove me?
71 validParamList->set< int > ("visualization: output file: iter", 0, "nonlinear iteration variable for output file name");//Remove me?
72 validParamList->set<std::string> ("visualization: style", "Point Cloud", "style of aggregate visualization for VTK output. Can be 'Point Cloud', 'Jacks', 'Convex Hulls'");
73 validParamList->set<bool> ("visualization: build colormap", false, "Whether to build a random color map in a separate xml file.");
74 validParamList->set<bool> ("visualization: fine graph edges", false, "Whether to draw all fine node connections along with the aggregates.");
75
76 return validParamList;
77 }
78
79 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::doPointCloud(std::vector<int>& vertices, std::vector<int>& geomSizes, LO /* numLocalAggs */, LO numFineNodes) {
81 vertices.reserve(numFineNodes);
82 geomSizes.reserve(numFineNodes);
83 for(LocalOrdinal i = 0; i < numFineNodes; i++)
84 {
85 vertices.push_back(i);
86 geomSizes.push_back(1);
87 }
88 }
89
90 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 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) {
92 //For each aggregate, find the root node then connect it with the other nodes in the aggregate
93 //Doesn't matter the order, as long as all edges are found.
94 vertices.reserve(vertices.size() + 3 * (numFineNodes - numLocalAggs));
95 geomSizes.reserve(vertices.size() + 2 * (numFineNodes - numLocalAggs));
96 int root = 0;
97 for(int i = 0; i < numLocalAggs; i++) //TODO: Replace this O(n^2) with a better way
98 {
99 while(!isRoot[root])
100 root++;
101 int numInAggFound = 0;
102 for(int j = 0; j < numFineNodes; j++)
103 {
104 if(j == root) //don't make a connection from the root to itself
105 {
106 numInAggFound++;
107 continue;
108 }
109 if(vertex2AggId[root] == vertex2AggId[j])
110 {
111 vertices.push_back(root);
112 vertices.push_back(j);
113 geomSizes.push_back(2);
114 //Also draw the free endpoint explicitly for the current line
115 vertices.push_back(j);
116 geomSizes.push_back(1);
117 numInAggFound++;
118 //if(numInAggFound == aggSizes_[vertex2AggId_[root]]) //don't spend more time looking if done with that root
119 // break;
120 }
121 }
122 root++; //get set up to look for the next root
123 }
124 }
125
126 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
127 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) {
129 // This algorithm is based on Andrew's Monotone Chain variant of the Graham Scan for Convex Hulls. It adds
130 // a colinearity check which we'll need for our viz.
131 for(int agg = 0; agg < numLocalAggs; agg++) {
132 std::vector<int> aggNodes;
133 for(int i = 0; i < numFineNodes; i++) {
134 if(vertex2AggId[i] == agg)
135 aggNodes.push_back(i);
138 //have a list of nodes in the aggregate
139 TEUCHOS_TEST_FOR_EXCEPTION(aggNodes.size() == 0, Exceptions::RuntimeError,
140 "CoarseningVisualization::doConvexHulls2D: aggregate contains zero nodes!");
141 if(aggNodes.size() == 1) {
142 vertices.push_back(aggNodes.front());
143 geomSizes.push_back(1);
144 continue;
145 }
146 if(aggNodes.size() == 2) {
147 vertices.push_back(aggNodes.front());
148 vertices.push_back(aggNodes.back());
149 geomSizes.push_back(2);
150 continue;
152 else {
153 int N = (int) aggNodes.size();
154 using MyPair = std::pair<myVec2,int>;
155 std::vector<MyPair> pointsAndIndex(N);
156 for(int i=0; i<N; i++) {
157 pointsAndIndex[i] = std::make_pair(myVec2(xCoords[aggNodes[i]], yCoords[aggNodes[i]]),i);
159
160 // Sort by x coordinate
161 std::sort(pointsAndIndex.begin(),pointsAndIndex.end(),[](const MyPair &a, const MyPair &b) {
162 return a.first.x < b.first.x || (a.first.x == b.first.x && a.first.y < b.first.y);
163 });
165
166 // Colinearity check
167 bool colinear=true;
168 for(int i=0; i<N; i++) {
169 if(ccw(pointsAndIndex[i].first,pointsAndIndex[(i+1)%N].first,pointsAndIndex[(i+2)%N].first)!=0) {
170 colinear=false;
171 break;
172 }
173 }
174
175 if(colinear) {
176 vertices.push_back(aggNodes[pointsAndIndex.front().second]);
177 vertices.push_back(aggNodes[pointsAndIndex.back().second]);
178 geomSizes.push_back(2);
179 }
180 else {
181 std::vector<int> hull(2*N);
182 int count=0;
183
184 // Build lower hull
185 for(int i=0; i<N; i++) {
186 while (count >=2 && ccw(pointsAndIndex[hull[count-2]].first,pointsAndIndex[hull[count-1]].first,pointsAndIndex[i].first)<=0) {
187 count--;
188 }
189 hull[count++] = i;
190 }
191
192 // Build the upper hull
193 int t=count+1;
194 for(int i=N-1; i>0; i--) {
195 while(count >= t && ccw(pointsAndIndex[hull[count-2]].first,pointsAndIndex[hull[count-1]].first,pointsAndIndex[i-1].first)<=0) {
196 count--;
197 }
198 hull[count++]=i-1;
199 }
200
201 // Remove the duplicated point
202 hull.resize(count-1);
203
204 // Verify: Check that hull retains CCW order
205 for(int i=0; i<(int)hull.size(); i++) {
206 TEUCHOS_TEST_FOR_EXCEPTION(ccw(pointsAndIndex[hull[i]].first,pointsAndIndex[hull[(i+1)%hull.size()]].first,pointsAndIndex[hull[(i+1)%hull.size()]].first) == 1, Exceptions::RuntimeError,"CoarseningVisualization::doConvexHulls2D: counter-clockwise verification fails");
207 }
208
209 // We now need to undo the indices from the sorting
210 for(int i=0; i<(int)hull.size(); i++) {
211 hull[i] = aggNodes[pointsAndIndex[hull[i]].second];
212 }
213
214 vertices.reserve(vertices.size() + hull.size());
215 for(size_t i = 0; i < hull.size(); i++) {
216 vertices.push_back(hull[i]);
217 }
218 geomSizes.push_back(hull.size());
219 }//else colinear
220 }//else 3 + nodes
221 }
222 }
223
224 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
225 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) {
226 //Use 3D quickhull algo.
227 //Vector of node indices representing triangle vertices
228 //Note: Calculate the hulls first since will only include point data for points in the hulls
229 //Effectively the size() of vertIndices after each hull is added to it
230 typedef std::list<int>::iterator Iter;
231 for(int agg = 0; agg < numLocalAggs; agg++) {
232 std::list<int> aggNodes; //At first, list of all nodes in the aggregate. As nodes are enclosed or included by/in hull, remove them
233 for(int i = 0; i < numFineNodes; i++) {
234 if(vertex2AggId[i] == agg)
235 aggNodes.push_back(i);
236 }
237 //First, check anomalous cases
238 TEUCHOS_TEST_FOR_EXCEPTION(aggNodes.size() == 0, Exceptions::RuntimeError,
239 "CoarseningVisualization::doConvexHulls3D: aggregate contains zero nodes!");
240 if(aggNodes.size() == 1) {
241 vertices.push_back(aggNodes.front());
242 geomSizes.push_back(1);
243 continue;
244 } else if(aggNodes.size() == 2) {
245 vertices.push_back(aggNodes.front());
246 vertices.push_back(aggNodes.back());
247 geomSizes.push_back(2);
248 continue;
249 }
250 //check for collinearity
251 bool areCollinear = true;
252 {
253 Iter it = aggNodes.begin();
254 myVec3 firstVec(xCoords[*it], yCoords[*it], zCoords[*it]);
255 myVec3 comp;
256 {
257 it++;
258 myVec3 secondVec(xCoords[*it], yCoords[*it], zCoords[*it]); //cross this with other vectors to compare
259 comp = vecSubtract(secondVec, firstVec);
260 it++;
261 }
262 for(; it != aggNodes.end(); it++) {
263 myVec3 thisVec(xCoords[*it], yCoords[*it], zCoords[*it]);
264 myVec3 cross = crossProduct(vecSubtract(thisVec, firstVec), comp);
265 if(mymagnitude(cross) > 1e-10) {
266 areCollinear = false;
267 break;
268 }
269 }
270 }
271 if(areCollinear)
272 {
273 //find the endpoints of segment describing all the points
274 //compare x, if tie compare y, if tie compare z
275 Iter min = aggNodes.begin();
276 Iter max = aggNodes.begin();
277 Iter it = ++aggNodes.begin();
278 for(; it != aggNodes.end(); it++) {
279 if(xCoords[*it] < xCoords[*min]) min = it;
280 else if(xCoords[*it] == xCoords[*min]) {
281 if(yCoords[*it] < yCoords[*min]) min = it;
282 else if(yCoords[*it] == yCoords[*min]) {
283 if(zCoords[*it] < zCoords[*min]) min = it;
284 }
285 }
286 if(xCoords[*it] > xCoords[*max]) max = it;
287 else if(xCoords[*it] == xCoords[*max]) {
288 if(yCoords[*it] > yCoords[*max]) max = it;
289 else if(yCoords[*it] == yCoords[*max]) {
290 if(zCoords[*it] > zCoords[*max])
291 max = it;
292 }
293 }
294 }
295 vertices.push_back(*min);
296 vertices.push_back(*max);
297 geomSizes.push_back(2);
298 continue;
299 }
300 bool areCoplanar = true;
301 {
302 //number of points is known to be >= 3
303 Iter vert = aggNodes.begin();
304 myVec3 v1(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
305 vert++;
306 myVec3 v2(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
307 vert++;
308 myVec3 v3(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
309 vert++;
310 //Make sure the first three points aren't also collinear (need a non-degenerate triangle to get a normal)
311 while(mymagnitude(crossProduct(vecSubtract(v1, v2), vecSubtract(v1, v3))) < 1e-10) {
312 //Replace the third point with the next point
313 v3 = myVec3(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
314 vert++;
315 }
316 for(; vert != aggNodes.end(); vert++) {
317 myVec3 pt(xCoords[*vert], yCoords[*vert], zCoords[*vert]);
318 if(fabs(pointDistFromTri(pt, v1, v2, v3)) > 1e-12) {
319 areCoplanar = false;
320 break;
321 }
322 }
323 if(areCoplanar) {
324 //do 2D convex hull
325 myVec3 planeNorm = getNorm(v1, v2, v3);
326 planeNorm.x = fabs(planeNorm.x);
327 planeNorm.y = fabs(planeNorm.y);
328 planeNorm.z = fabs(planeNorm.z);
329 std::vector<myVec2> points;
330 std::vector<int> nodes;
331 if(planeNorm.x >= planeNorm.y && planeNorm.x >= planeNorm.z) {
332 //project points to yz plane to make hull
333 for(Iter it = aggNodes.begin(); it != aggNodes.end(); it++) {
334 nodes.push_back(*it);
335 points.push_back(myVec2(yCoords[*it], zCoords[*it]));
336 }
337 }
338 if(planeNorm.y >= planeNorm.x && planeNorm.y >= planeNorm.z) {
339 //use xz
340 for(Iter it = aggNodes.begin(); it != aggNodes.end(); it++) {
341 nodes.push_back(*it);
342 points.push_back(myVec2(xCoords[*it], zCoords[*it]));
343 }
344 }
345 if(planeNorm.z >= planeNorm.x && planeNorm.z >= planeNorm.y) {
346 for(Iter it = aggNodes.begin(); it != aggNodes.end(); it++) {
347 nodes.push_back(*it);
348 points.push_back(myVec2(xCoords[*it], yCoords[*it]));
349 }
350 }
351 std::vector<int> convhull2d = giftWrap(points, nodes, xCoords, yCoords);
352 geomSizes.push_back(convhull2d.size());
353 vertices.reserve(vertices.size() + convhull2d.size());
354 for(size_t i = 0; i < convhull2d.size(); i++)
355 vertices.push_back(convhull2d[i]);
356 continue;
357 }
358 }
359 Iter exIt = aggNodes.begin(); //iterator to be used for searching for min/max x/y/z
360 int extremeSix[] = {*exIt, *exIt, *exIt, *exIt, *exIt, *exIt}; //nodes with minimumX, maxX, minY ...
361 exIt++;
362 for(; exIt != aggNodes.end(); exIt++) {
363 Iter it = exIt;
364 if(xCoords[*it] < xCoords[extremeSix[0]] ||
365 (xCoords[*it] == xCoords[extremeSix[0]] && yCoords[*it] < yCoords[extremeSix[0]]) ||
366 (xCoords[*it] == xCoords[extremeSix[0]] && yCoords[*it] == yCoords[extremeSix[0]] && zCoords[*it] < zCoords[extremeSix[0]]))
367 extremeSix[0] = *it;
368 if(xCoords[*it] > xCoords[extremeSix[1]] ||
369 (xCoords[*it] == xCoords[extremeSix[1]] && yCoords[*it] > yCoords[extremeSix[1]]) ||
370 (xCoords[*it] == xCoords[extremeSix[1]] && yCoords[*it] == yCoords[extremeSix[1]] && zCoords[*it] > zCoords[extremeSix[1]]))
371 extremeSix[1] = *it;
372 if(yCoords[*it] < yCoords[extremeSix[2]] ||
373 (yCoords[*it] == yCoords[extremeSix[2]] && zCoords[*it] < zCoords[extremeSix[2]]) ||
374 (yCoords[*it] == yCoords[extremeSix[2]] && zCoords[*it] == zCoords[extremeSix[2]] && xCoords[*it] < xCoords[extremeSix[2]]))
375 extremeSix[2] = *it;
376 if(yCoords[*it] > yCoords[extremeSix[3]] ||
377 (yCoords[*it] == yCoords[extremeSix[3]] && zCoords[*it] > zCoords[extremeSix[3]]) ||
378 (yCoords[*it] == yCoords[extremeSix[3]] && zCoords[*it] == zCoords[extremeSix[3]] && xCoords[*it] > xCoords[extremeSix[3]]))
379 extremeSix[3] = *it;
380 if(zCoords[*it] < zCoords[extremeSix[4]] ||
381 (zCoords[*it] == zCoords[extremeSix[4]] && xCoords[*it] < xCoords[extremeSix[4]]) ||
382 (zCoords[*it] == zCoords[extremeSix[4]] && xCoords[*it] == xCoords[extremeSix[4]] && yCoords[*it] < yCoords[extremeSix[4]]))
383 extremeSix[4] = *it;
384 if(zCoords[*it] > zCoords[extremeSix[5]] ||
385 (zCoords[*it] == zCoords[extremeSix[5]] && xCoords[*it] > xCoords[extremeSix[5]]) ||
386 (zCoords[*it] == zCoords[extremeSix[5]] && xCoords[*it] == xCoords[extremeSix[5]] && yCoords[*it] > zCoords[extremeSix[5]]))
387 extremeSix[5] = *it;
388 }
389 myVec3 extremeVectors[6];
390 for(int i = 0; i < 6; i++) {
391 myVec3 thisExtremeVec(xCoords[extremeSix[i]], yCoords[extremeSix[i]], zCoords[extremeSix[i]]);
392 extremeVectors[i] = thisExtremeVec;
393 }
394 double maxDist = 0;
395 int base1 = 0; //ints from 0-5: which pair out of the 6 extreme points are the most distant? (indices in extremeSix and extremeVectors)
396 int base2 = 0;
397 for(int i = 0; i < 5; i++) {
398 for(int j = i + 1; j < 6; j++) {
399 double thisDist = distance(extremeVectors[i], extremeVectors[j]);
400 // Want to make sure thisDist is significantly larger than maxDist to prevent roundoff errors from impacting node choice.
401 if(thisDist > maxDist && thisDist - maxDist > 1e-12) {
402 maxDist = thisDist;
403 base1 = i;
404 base2 = j;
405 }
406 }
407 }
408 std::list<myTriangle> hullBuilding; //each Triangle is a triplet of nodes (int IDs) that form a triangle
409 //remove base1 and base2 iters from aggNodes, they are known to be in the hull
410 aggNodes.remove(extremeSix[base1]);
411 aggNodes.remove(extremeSix[base2]);
412 //extremeSix[base1] and [base2] still have the myVec3 representation
413 myTriangle tri;
414 tri.v1 = extremeSix[base1];
415 tri.v2 = extremeSix[base2];
416 //Now find the point that is furthest away from the line between base1 and base2
417 maxDist = 0;
418 //need the vectors to do "quadruple product" formula
419 myVec3 b1 = extremeVectors[base1];
420 myVec3 b2 = extremeVectors[base2];
421 Iter thirdNode;
422 for(Iter node = aggNodes.begin(); node != aggNodes.end(); node++) {
423 myVec3 nodePos(xCoords[*node], yCoords[*node], zCoords[*node]);
424 double dist = mymagnitude(crossProduct(vecSubtract(nodePos, b1), vecSubtract(nodePos, b2))) / mymagnitude(vecSubtract(b2, b1));
425 // Want to make sure dist is significantly larger than maxDist to prevent roundoff errors from impacting node choice.
426 if(dist > maxDist && dist - maxDist > 1e-12) {
427 maxDist = dist;
428 thirdNode = node;
429 }
430 }
431 //Now know the last node in the first triangle
432 tri.v3 = *thirdNode;
433 hullBuilding.push_back(tri);
434 myVec3 b3(xCoords[*thirdNode], yCoords[*thirdNode], zCoords[*thirdNode]);
435 aggNodes.erase(thirdNode);
436 //Find the fourth node (most distant from triangle) to form tetrahedron
437 maxDist = 0;
438 int fourthVertex = -1;
439 for(Iter node = aggNodes.begin(); node != aggNodes.end(); node++) {
440 myVec3 thisNode(xCoords[*node], yCoords[*node], zCoords[*node]);
441 double nodeDist = pointDistFromTri(thisNode, b1, b2, b3);
442 // Want to make sure nodeDist is significantly larger than maxDist to prevent roundoff errors from impacting node choice.
443 if(nodeDist > maxDist && nodeDist - maxDist > 1e-12) {
444 maxDist = nodeDist;
445 fourthVertex = *node;
446 }
447 }
448 aggNodes.remove(fourthVertex);
449 myVec3 b4(xCoords[fourthVertex], yCoords[fourthVertex], zCoords[fourthVertex]);
450 //Add three new triangles to hullBuilding to form the first tetrahedron
451 //use tri to hold the triangle info temporarily before being added to list
452 tri = hullBuilding.front();
453 tri.v1 = fourthVertex;
454 hullBuilding.push_back(tri);
455 tri = hullBuilding.front();
456 tri.v2 = fourthVertex;
457 hullBuilding.push_back(tri);
458 tri = hullBuilding.front();
459 tri.v3 = fourthVertex;
460 hullBuilding.push_back(tri);
461 //now orient all four triangles so that the vertices are oriented clockwise (so getNorm_ points outward for each)
462 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);
463 for(std::list<myTriangle>::iterator tetTri = hullBuilding.begin(); tetTri != hullBuilding.end(); tetTri++) {
464 myVec3 triVert1(xCoords[tetTri->v1], yCoords[tetTri->v1], zCoords[tetTri->v1]);
465 myVec3 triVert2(xCoords[tetTri->v2], yCoords[tetTri->v2], zCoords[tetTri->v2]);
466 myVec3 triVert3(xCoords[tetTri->v3], yCoords[tetTri->v3], zCoords[tetTri->v3]);
467 myVec3 trinorm = getNorm(triVert1, triVert2, triVert3);
468 myVec3 ptInPlane = tetTri == hullBuilding.begin() ? b1 : b4; //first triangle definitely has b1 in it, other three definitely have b4
469 if(isInFront(barycenter, ptInPlane, trinorm)) {
470 //don't want the faces of the tetrahedron to face inwards (towards barycenter) so reverse orientation
471 //by swapping two vertices
472 int temp = tetTri->v1;
473 tetTri->v1 = tetTri->v2;
474 tetTri->v2 = temp;
475 }
476 }
477 //now, have starting polyhedron in hullBuilding (all faces are facing outwards according to getNorm_) and remaining nodes to process are in aggNodes
478 //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.
479 //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
480 //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
481 //new triangles, since they are now enclosed in the hull.
482 //Construct point lists for each face of the tetrahedron individually.
483 myVec3 trinorms[4]; //normals to the four tetrahedron faces, now oriented outwards
484 int index = 0;
485 for(std::list<myTriangle>::iterator tetTri = hullBuilding.begin(); tetTri != hullBuilding.end(); tetTri++) {
486 myVec3 triVert1(xCoords[tetTri->v1], yCoords[tetTri->v1], zCoords[tetTri->v1]);
487 myVec3 triVert2(xCoords[tetTri->v2], yCoords[tetTri->v2], zCoords[tetTri->v2]);
488 myVec3 triVert3(xCoords[tetTri->v3], yCoords[tetTri->v3], zCoords[tetTri->v3]);
489 trinorms[index] = getNorm(triVert1, triVert2, triVert3);
490 index++;
491 }
492 std::list<int> startPoints1;
493 std::list<int> startPoints2;
494 std::list<int> startPoints3;
495 std::list<int> startPoints4;
496 //scope this so that 'point' is not in function scope
497 {
498 Iter point = aggNodes.begin();
499 while(!aggNodes.empty()) //this removes points one at a time as they are put in startPointsN or are already done
500 {
501 myVec3 pointVec(xCoords[*point], yCoords[*point], zCoords[*point]);
502 //Note: Because of the way the tetrahedron faces are constructed above,
503 //face 1 definitely contains b1 and faces 2-4 definitely contain b4.
504 if(isInFront(pointVec, b1, trinorms[0])) {
505 startPoints1.push_back(*point);
506 point = aggNodes.erase(point);
507 } else if(isInFront(pointVec, b4, trinorms[1])) {
508 startPoints2.push_back(*point);
509 point = aggNodes.erase(point);
510 } else if(isInFront(pointVec, b4, trinorms[2])) {
511 startPoints3.push_back(*point);
512 point = aggNodes.erase(point);
513 } else if(isInFront(pointVec, b4, trinorms[3])) {
514 startPoints4.push_back(*point);
515 point = aggNodes.erase(point);
516 } else {
517 point = aggNodes.erase(point); //points here are already inside tetrahedron.
518 }
519 }
520 //Call processTriangle for each triangle in the initial tetrahedron, one at a time.
521 }
522 typedef std::list<myTriangle>::iterator TriIter;
523 TriIter firstTri = hullBuilding.begin();
524 myTriangle start1 = *firstTri;
525 firstTri++;
526 myTriangle start2 = *firstTri;
527 firstTri++;
528 myTriangle start3 = *firstTri;
529 firstTri++;
530 myTriangle start4 = *firstTri;
531 //kick off depth-first recursive filling of hullBuilding list with all triangles in the convex hull
532 if(!startPoints1.empty())
533 processTriangle(hullBuilding, start1, startPoints1, barycenter, xCoords, yCoords, zCoords);
534 if(!startPoints2.empty())
535 processTriangle(hullBuilding, start2, startPoints2, barycenter, xCoords, yCoords, zCoords);
536 if(!startPoints3.empty())
537 processTriangle(hullBuilding, start3, startPoints3, barycenter, xCoords, yCoords, zCoords);
538 if(!startPoints4.empty())
539 processTriangle(hullBuilding, start4, startPoints4, barycenter, xCoords, yCoords, zCoords);
540 //hullBuilding now has all triangles that make up this hull.
541 //Dump hullBuilding info into the list of all triangles for the scene.
542 vertices.reserve(vertices.size() + 3 * hullBuilding.size());
543 for(TriIter hullTri = hullBuilding.begin(); hullTri != hullBuilding.end(); hullTri++) {
544 vertices.push_back(hullTri->v1);
545 vertices.push_back(hullTri->v2);
546 vertices.push_back(hullTri->v3);
547 geomSizes.push_back(3);
548 }
549 }
550 }
551
552 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
553 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::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 */) {
554 ArrayView<const Scalar> values;
555 ArrayView<const LocalOrdinal> neighbors;
556
557 std::vector<std::pair<int, int> > vert1; //vertices (node indices)
558
559 ArrayView<const LocalOrdinal> indices;
560 for(LocalOrdinal locRow = 0; locRow < LocalOrdinal(G->GetNodeNumVertices()); locRow++) {
561 neighbors = G->getNeighborVertices(locRow);
562 //Add those local indices (columns) to the list of connections (which are pairs of the form (localM, localN))
563 for(int gEdge = 0; gEdge < int(neighbors.size()); ++gEdge) {
564 vert1.push_back(std::pair<int, int>(locRow, neighbors[gEdge]));
565 }
566 }
567 for(size_t i = 0; i < vert1.size(); i ++) {
568 if(vert1[i].first > vert1[i].second) {
569 int temp = vert1[i].first;
570 vert1[i].first = vert1[i].second;
571 vert1[i].second = temp;
572 }
573 }
574 std::sort(vert1.begin(), vert1.end());
575 std::vector<std::pair<int, int> >::iterator newEnd = unique(vert1.begin(), vert1.end()); //remove duplicate edges
576 vert1.erase(newEnd, vert1.end());
577 //std::vector<int> points1;
578 vertices.reserve(2 * vert1.size());
579 geomSizes.reserve(vert1.size());
580 for(size_t i = 0; i < vert1.size(); i++) {
581 vertices.push_back(vert1[i].first);
582 vertices.push_back(vert1[i].second);
583 geomSizes.push_back(2);
584 }
585 }
586
587 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
589 {
590 const double ccw_zero_threshold=1e-8;
591 double val = (b.x - a.x) * (c.y - a.y) - (b.y - a.y) * (c.x - a.x);
592 return (val > ccw_zero_threshold) ? 1 : ( (val < -ccw_zero_threshold) ? -1 : 0);
593 }
594
595 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
597 {
598 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);
599 }
600
601 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
603 {
604 return v1.x * v2.x + v1.y * v2.y;
605 }
606
607 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
609 {
610 return v1.x * v2.x + v1.y * v2.y + v1.z * v2.z;
611 }
612
613 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
615 {
616 myVec3 rel(point.x - inPlane.x, point.y - inPlane.y, point.z - inPlane.z); //position of the point relative to the plane
617 return dotProduct(rel, n) > 1e-12 ? true : false;
618 }
619
620 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
625
626 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
631
632 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
637
638
639 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
644
645 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
650
651 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
656
657 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
658 myVec2 VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getNorm(myVec2 v) //"normal" to a 2D vector - just rotate 90 degrees to left
659 {
660 return myVec2(v.y, -v.x);
661 }
662
663 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
664 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)
665 {
666 return crossProduct(vecSubtract(v2, v1), vecSubtract(v3, v1));
667 }
668
669 //get minimum distance from 'point' to plane containing v1, v2, v3 (or the triangle with v1, v2, v3 as vertices)
670 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
672 {
673 using namespace std;
674 myVec3 norm = getNorm(v1, v2, v3);
675 //must normalize the normal vector
676 double normScl = mymagnitude(norm);
677 double rv = 0.0;
678 if (normScl > 1e-8) {
679 norm.x /= normScl;
680 norm.y /= normScl;
681 norm.z /= normScl;
682 rv = fabs(dotProduct(norm, vecSubtract(point, v1)));
683 } else {
684 // triangle is degenerated
685 myVec3 test1 = vecSubtract(v3, v1);
686 myVec3 test2 = vecSubtract(v2, v1);
687 bool useTest1 = mymagnitude(test1) > 0.0 ? true : false;
688 bool useTest2 = mymagnitude(test2) > 0.0 ? true : false;
689 if(useTest1 == true) {
690 double normScl1 = mymagnitude(test1);
691 test1.x /= normScl1;
692 test1.y /= normScl1;
693 test1.z /= normScl1;
694 rv = fabs(dotProduct(test1, vecSubtract(point, v1)));
695 } else if (useTest2 == true) {
696 double normScl2 = mymagnitude(test2);
697 test2.x /= normScl2;
698 test2.y /= normScl2;
699 test2.z /= normScl2;
700 rv = fabs(dotProduct(test2, vecSubtract(point, v1)));
701 } else {
702 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError,
703 "VisualizationHelpers::pointDistFromTri: Could not determine the distance of a point to a triangle.");
704 }
705
706 }
707 return rv;
708 }
709
710 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
711 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) {
712 //*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.
713 //precondition: each triangle is already oriented so that getNorm_(v1, v2, v3) points outward (away from interior of hull)
714 //First find the point furthest from triangle.
715 using namespace std;
716 typedef std::list<int>::iterator Iter;
717 typedef std::list<myTriangle>::iterator TriIter;
718 typedef list<pair<int, int> >::iterator EdgeIter;
719 double maxDist = 0;
720 //Need vector representations of triangle's vertices
721 myVec3 v1(xCoords[tri.v1], yCoords[tri.v1], zCoords[tri.v1]);
722 myVec3 v2(xCoords[tri.v2], yCoords[tri.v2], zCoords[tri.v2]);
723 myVec3 v3(xCoords[tri.v3], yCoords[tri.v3], zCoords[tri.v3]);
724 myVec3 farPointVec; //useful to have both the point's coordinates and it's position in the list
725 Iter farPoint = pointsInFront.begin();
726 for(Iter point = pointsInFront.begin(); point != pointsInFront.end(); point++)
727 {
728 myVec3 pointVec(xCoords[*point], yCoords[*point], zCoords[*point]);
729 double dist = pointDistFromTri(pointVec, v1, v2, v3);
730 // Want to make sure nodeDist is significantly larger than maxDist to prevent roundoff errors from impacting node choice.
731 if(dist > maxDist && dist - maxDist > 1e-12)
732 {
733 dist = maxDist;
734 farPointVec = pointVec;
735 farPoint = point;
736 }
737 }
738 //Find all the triangles that the point is in front of (can be more than 1)
739 //At the same time, remove them from tris, as every one will be replaced later
740 vector<myTriangle> visible; //use a list of iterators so that the underlying object is still in tris
741 for(TriIter it = tris.begin(); it != tris.end();)
742 {
743 myVec3 vec1(xCoords[it->v1], yCoords[it->v1], zCoords[it->v1]);
744 myVec3 vec2(xCoords[it->v2], yCoords[it->v2], zCoords[it->v2]);
745 myVec3 vec3(xCoords[it->v3], yCoords[it->v3], zCoords[it->v3]);
746 myVec3 norm = getNorm(vec1, vec2, vec3);
747 if(isInFront(farPointVec, vec1, norm))
748 {
749 visible.push_back(*it);
750 it = tris.erase(it);
751 }
752 else
753 it++;
754 }
755 //Figure out what triangles need to be destroyed/created
756 //First create a list of edges (as std::pair<int, int>, where the two ints are the node endpoints)
757 list<pair<int, int> > horizon;
758 //For each triangle, add edges to the list iff the edge only appears once in the set of all
759 //Have members of horizon have the lower node # first, and the higher one second
760 for(vector<myTriangle>::iterator it = visible.begin(); it != visible.end(); it++)
761 {
762 pair<int, int> e1(it->v1, it->v2);
763 pair<int, int> e2(it->v2, it->v3);
764 pair<int, int> e3(it->v1, it->v3);
765 //"sort" the pair values
766 if(e1.first > e1.second)
767 {
768 int temp = e1.first;
769 e1.first = e1.second;
770 e1.second = temp;
771 }
772 if(e2.first > e2.second)
773 {
774 int temp = e2.first;
775 e2.first = e2.second;
776 e2.second = temp;
777 }
778 if(e3.first > e3.second)
779 {
780 int temp = e3.first;
781 e3.first = e3.second;
782 e3.second = temp;
783 }
784 horizon.push_back(e1);
785 horizon.push_back(e2);
786 horizon.push_back(e3);
787 }
788 //sort based on lower node first, then higher node (lexicographical ordering built in to pair)
789 horizon.sort();
790 //Remove all edges from horizon, except those that appear exactly once
791 {
792 EdgeIter it = horizon.begin();
793 while(it != horizon.end())
794 {
795 int occur = count(horizon.begin(), horizon.end(), *it);
796 if(occur > 1)
797 {
798 pair<int, int> removeVal = *it;
799 while(removeVal == *it && !(it == horizon.end()))
800 it = horizon.erase(it);
801 }
802 else
803 it++;
804 }
805 }
806 //Now make a list of new triangles being created, each of which take 2 vertices from an edge and one from farPoint
807 list<myTriangle> newTris;
808 for(EdgeIter it = horizon.begin(); it != horizon.end(); it++)
809 {
810 myTriangle t(it->first, it->second, *farPoint);
811 newTris.push_back(t);
812 }
813 //Ensure every new triangle is oriented outwards, using the barycenter of the initial tetrahedron
814 vector<myTriangle> trisToProcess;
815 vector<list<int> > newFrontPoints;
816 for(TriIter it = newTris.begin(); it != newTris.end(); it++)
817 {
818 myVec3 t1(xCoords[it->v1], yCoords[it->v1], zCoords[it->v1]);
819 myVec3 t2(xCoords[it->v2], yCoords[it->v2], zCoords[it->v2]);
820 myVec3 t3(xCoords[it->v3], yCoords[it->v3], zCoords[it->v3]);
821 if(isInFront(barycenter, t1, getNorm(t1, t2, t3)))
822 {
823 //need to swap two vertices to flip orientation of triangle
824 int temp = it->v1;
825 myVec3 tempVec = t1;
826 it->v1 = it->v2;
827 t1 = t2;
828 it->v2 = temp;
829 t2 = tempVec;
830 }
831 myVec3 outwardNorm = getNorm(t1, t2, t3); //now definitely points outwards
832 //Add the triangle to tris
833 tris.push_back(*it);
834 trisToProcess.push_back(tris.back());
835 //Make a list of the points that are in front of nextToProcess, to be passed in for processing
836 list<int> newInFront;
837 for(Iter point = pointsInFront.begin(); point != pointsInFront.end();)
838 {
839 myVec3 pointVec(xCoords[*point], yCoords[*point], zCoords[*point]);
840 if(isInFront(pointVec, t1, outwardNorm))
841 {
842 newInFront.push_back(*point);
843 point = pointsInFront.erase(point);
844 }
845 else
846 point++;
847 }
848 newFrontPoints.push_back(newInFront);
849 }
850 vector<myTriangle> allRemoved; //list of all invalid iterators that were erased by calls to processmyTriangle below
851 for(int i = 0; i < int(trisToProcess.size()); i++)
852 {
853 if(!newFrontPoints[i].empty())
854 {
855 //Comparing the 'triangle to process' to the one for this call prevents infinite recursion/stack overflow.
856 //TODO: Why was it doing that? Rounding error? Make more robust fix. But this does work for the time being.
857 if(find(allRemoved.begin(), allRemoved.end(), trisToProcess[i]) == allRemoved.end() && !(trisToProcess[i] == tri))
858 {
859 vector<myTriangle> removedList = processTriangle(tris, trisToProcess[i], newFrontPoints[i], barycenter, xCoords, yCoords, zCoords);
860 for(int j = 0; j < int(removedList.size()); j++)
861 allRemoved.push_back(removedList[j]);
862 }
863 }
864 }
865 return visible;
866 }
867
868 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
869 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) {
870 TEUCHOS_TEST_FOR_EXCEPTION(points.size() < 3, Exceptions::RuntimeError,
871 "CoarseningVisualization::giftWrap: Gift wrap algorithm input has to have at least 3 points!");
872
873#if 1 // TAW's version to determine "minimal" node
874 // determine minimal x and y coordinates
875 double min_x =points[0].x;
876 double min_y =points[0].y;
877 for(std::vector<int>::iterator it = nodes.begin(); it != nodes.end(); it++) {
878 int i = it - nodes.begin();
879 if(points[i].x < min_x) min_x = points[i].x;
880 if(points[i].y < min_y) min_y = points[i].y;
881 }
882 // create dummy min coordinates
883 min_x -= 1.0;
884 min_y -= 1.0;
885 myVec2 dummy_min(min_x, min_y);
886
887 // loop over all nodes and determine nodes with minimal distance to (min_x, min_y)
888 std::vector<int> hull;
889 myVec2 min = points[0];
890 double mindist = distance(min,dummy_min);
891 std::vector<int>::iterator minNode = nodes.begin();
892 for(std::vector<int>::iterator it = nodes.begin(); it != nodes.end(); it++) {
893 int i = it - nodes.begin();
894 if(distance(points[i],dummy_min) < mindist) {
895 mindist = distance(points[i],dummy_min);
896 min = points[i];
897 minNode = it;
898 }
899 }
900 hull.push_back(*minNode);
901#else // Brian's code
902 std::vector<int> hull;
903 std::vector<int>::iterator minNode = nodes.begin();
904 myVec2 min = points[0];
905 for(std::vector<int>::iterator it = nodes.begin(); it != nodes.end(); it++)
906 {
907 int i = it - nodes.begin();
908 if(points[i].x < min.x || (fabs(points[i].x - min.x) < 1e-10 && points[i].y < min.y))
909 {
910 min = points[i];
911 minNode = it;
912 }
913 }
914 hull.push_back(*minNode);
915#endif
916
917 bool includeMin = false;
918 //int debug_it = 0;
919 while(1)
920 {
921 std::vector<int>::iterator leftMost = nodes.begin();
922 if(!includeMin && leftMost == minNode)
923 {
924 leftMost++;
925 }
926 std::vector<int>::iterator it = leftMost;
927 it++;
928 for(; it != nodes.end(); it++)
929 {
930 if(it == minNode && !includeMin) //don't compare to min on very first sweep
931 continue;
932 if(*it == hull.back())
933 continue;
934 //see if it is in front of line containing nodes thisHull.back() and leftMost
935 //first get the left normal of leftMost - thisHull.back() (<dy, -dx>)
936 myVec2 leftMostVec = points[leftMost - nodes.begin()];
937 myVec2 lastVec(xCoords[hull.back()], yCoords[hull.back()]);
938 myVec2 testNorm = getNorm(vecSubtract(leftMostVec, lastVec));
939 //now dot testNorm with *it - leftMost. If dot is positive, leftMost becomes it. If dot is zero, take one further from thisHull.back().
940 myVec2 itVec(xCoords[*it], yCoords[*it]);
941 double dotProd = dotProduct(testNorm, vecSubtract(itVec, lastVec));
942 if(-1e-8 < dotProd && dotProd < 1e-8)
943 {
944 //thisHull.back(), it and leftMost are collinear.
945 //Just sum the differences in x and differences in y for each and compare to get further one, don't need distance formula
946 myVec2 itDist = vecSubtract(itVec, lastVec);
947 myVec2 leftMostDist = vecSubtract(leftMostVec, lastVec);
948 if(fabs(itDist.x) + fabs(itDist.y) > fabs(leftMostDist.x) + fabs(leftMostDist.y)) {
949 leftMost = it;
950 }
951 }
952 else if(dotProd > 0) {
953 leftMost = it;
954
955 }
956 }
957 //if leftMost is min, then the loop is complete.
958 if(*leftMost == *minNode)
959 break;
960 hull.push_back(*leftMost);
961 includeMin = true; //have found second point (the one after min) so now include min in the searches
962 //debug_it ++;
963 //if(debug_it > 100) exit(0); //break;
964 }
965 return hull;
966 }
967
968 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
970 {
971 using namespace std;
972 vector<int> uniqueNodes = vertices;
973 sort(uniqueNodes.begin(), uniqueNodes.end());
974 vector<int>::iterator newUniqueFineEnd = unique(uniqueNodes.begin(), uniqueNodes.end());
975 uniqueNodes.erase(newUniqueFineEnd, uniqueNodes.end());
976 //uniqueNodes is now a sorted list of the nodes whose info actually goes in file
977 //Now replace values in vertices with locations of the old values in uniqueFine
978 for(int i = 0; i < int(vertices.size()); i++)
979 {
980 int lo = 0;
981 int hi = uniqueNodes.size() - 1;
982 int mid = 0;
983 int search = vertices[i];
984 while(lo <= hi)
985 {
986 mid = lo + (hi - lo) / 2;
987 if(uniqueNodes[mid] == search)
988 break;
989 else if(uniqueNodes[mid] > search)
990 hi = mid - 1;
991 else
992 lo = mid + 1;
993 }
994 if(uniqueNodes[mid] != search)
995 throw runtime_error("Issue in makeUnique_() - a point wasn't found in list.");
996 vertices[i] = mid;
997 }
998 return uniqueNodes;
999 }
1000
1001 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1002 std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceAll(std::string result, const std::string& replaceWhat, const std::string& replaceWithWhat) const {
1003 while(1) {
1004 const int pos = result.find(replaceWhat);
1005 if (pos == -1)
1006 break;
1007 result.replace(pos, replaceWhat.size(), replaceWithWhat);
1008 }
1009 return result;
1010 }
1011
1012 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1013 std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getFileName(int numProcs, int myRank, int level, const Teuchos::ParameterList & pL) const {
1014 std::string filenameToWrite = getBaseFileName(numProcs, level, pL);
1015 filenameToWrite = this->replaceAll(filenameToWrite, "%PROCID", toString(myRank));
1016 return filenameToWrite;
1017 }
1018
1019 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1020 std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getBaseFileName(int numProcs, int level, const Teuchos::ParameterList & pL) const {
1021 std::string filenameToWrite = pL.get<std::string>("visualization: output filename");
1022 int timeStep = pL.get<int>("visualization: output file: time step");
1023 int iter = pL.get<int>("visualization: output file: iter");
1024
1025 if(filenameToWrite.rfind(".vtu") == std::string::npos)
1026 filenameToWrite.append(".vtu");
1027 if(numProcs > 1 && filenameToWrite.rfind("%PROCID") == std::string::npos) //filename can't be identical between processsors in parallel problem
1028 filenameToWrite.insert(filenameToWrite.rfind(".vtu"), "-proc%PROCID");
1029
1030 filenameToWrite = this->replaceAll(filenameToWrite, "%LEVELID", toString(level));
1031 filenameToWrite = this->replaceAll(filenameToWrite, "%TIMESTEP", toString(timeStep));
1032 filenameToWrite = this->replaceAll(filenameToWrite, "%ITER", toString(iter));
1033 return filenameToWrite;
1034 }
1035
1036 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1037 std::string VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getPVTUFileName(int numProcs, int /* myRank */, int level, const Teuchos::ParameterList & pL) const {
1038 std::string filenameToWrite = getBaseFileName(numProcs, level, pL);
1039 std::string masterStem = filenameToWrite.substr(0, filenameToWrite.rfind(".vtu"));
1040 masterStem = this->replaceAll(masterStem, "%PROCID", "");
1041 std::string pvtuFilename = masterStem + "-master.pvtu";
1042 return pvtuFilename;
1043 }
1044
1045 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1046 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKOpening(std::ofstream & fout, std::vector<int> & uniqueFine, std::vector<int> & geomSizesFine) const {
1047 std::string styleName = "PointCloud"; // TODO adapt this
1048
1049 std::string indent = " ";
1050 fout << "<!--" << styleName << " Aggregates Visualization-->" << std::endl;
1051 fout << "<VTKFile type=\"UnstructuredGrid\" byte_order=\"LittleEndian\">" << std::endl;
1052 fout << " <UnstructuredGrid>" << std::endl;
1053 fout << " <Piece NumberOfPoints=\"" << uniqueFine.size() << "\" NumberOfCells=\"" << geomSizesFine.size() << "\">" << std::endl;
1054 fout << " <PointData Scalars=\"Node Aggregate Processor\">" << std::endl;
1055 }
1056
1057 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1058 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKNodes(std::ofstream & fout, std::vector<int> & uniqueFine, Teuchos::RCP<const Map> & nodeMap) const {
1059 std::string indent = " ";
1060 fout << " <DataArray type=\"Int32\" Name=\"Node\" format=\"ascii\">" << std::endl;
1061 indent = " ";
1062 bool localIsGlobal = GlobalOrdinal(nodeMap->getGlobalNumElements()) == GlobalOrdinal(nodeMap->getLocalNumElements());
1063 for(size_t i = 0; i < uniqueFine.size(); i++)
1064 {
1065 if(localIsGlobal)
1066 {
1067 fout << uniqueFine[i] << " "; //if all nodes are on this processor, do not map from local to global
1068 }
1069 else
1070 fout << nodeMap->getGlobalElement(uniqueFine[i]) << " ";
1071 if(i % 10 == 9)
1072 fout << std::endl << indent;
1073 }
1074 fout << std::endl;
1075 fout << " </DataArray>" << std::endl;
1076 }
1077
1078 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1079 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKData(std::ofstream & fout, std::vector<int> & uniqueFine, LocalOrdinal myAggOffset, ArrayRCP<LocalOrdinal> & vertex2AggId, int myRank) const {
1080 std::string indent = " ";
1081 fout << " <DataArray type=\"Int32\" Name=\"Aggregate\" format=\"ascii\">" << std::endl;
1082 fout << indent;
1083 for(int i = 0; i < int(uniqueFine.size()); i++)
1084 {
1085 fout << myAggOffset + vertex2AggId[uniqueFine[i]] << " ";
1086 if(i % 10 == 9)
1087 fout << std::endl << indent;
1088 }
1089 fout << std::endl;
1090 fout << " </DataArray>" << std::endl;
1091 fout << " <DataArray type=\"Int32\" Name=\"Processor\" format=\"ascii\">" << std::endl;
1092 fout << indent;
1093 for(int i = 0; i < int(uniqueFine.size()); i++)
1094 {
1095 fout << myRank << " ";
1096 if(i % 20 == 19)
1097 fout << std::endl << indent;
1098 }
1099 fout << std::endl;
1100 fout << " </DataArray>" << std::endl;
1101 fout << " </PointData>" << std::endl;
1102 }
1103
1104 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1105 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::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 {
1106 std::string indent = " ";
1107 fout << " <Points>" << std::endl;
1108 fout << " <DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">" << std::endl;
1109 fout << indent;
1110 for(int i = 0; i < int(uniqueFine.size()); i++)
1111 {
1112 fout << fx[uniqueFine[i]] << " " << fy[uniqueFine[i]] << " ";
1113 if(dim == 2)
1114 fout << "0 ";
1115 else
1116 fout << fz[uniqueFine[i]] << " ";
1117 if(i % 3 == 2)
1118 fout << std::endl << indent;
1119 }
1120 fout << std::endl;
1121 fout << " </DataArray>" << std::endl;
1122 fout << " </Points>" << std::endl;
1123 }
1124
1125 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1126 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writeFileVTKCells(std::ofstream & fout, std::vector<int> & /* uniqueFine */, std::vector<LocalOrdinal> & vertices, std::vector<LocalOrdinal> & geomSize) const {
1127 std::string indent = " ";
1128 fout << " <Cells>" << std::endl;
1129 fout << " <DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">" << std::endl;
1130 fout << indent;
1131 for(int i = 0; i < int(vertices.size()); i++)
1132 {
1133 fout << vertices[i] << " ";
1134 if(i % 10 == 9)
1135 fout << std::endl << indent;
1136 }
1137 fout << std::endl;
1138 fout << " </DataArray>" << std::endl;
1139 fout << " <DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">" << std::endl;
1140 fout << indent;
1141 int accum = 0;
1142 for(int i = 0; i < int(geomSize.size()); i++)
1143 {
1144 accum += geomSize[i];
1145 fout << accum << " ";
1146 if(i % 10 == 9)
1147 fout << std::endl << indent;
1148 }
1149 fout << std::endl;
1150 fout << " </DataArray>" << std::endl;
1151 fout << " <DataArray type=\"Int32\" Name=\"types\" format=\"ascii\">" << std::endl;
1152 fout << indent;
1153 for(int i = 0; i < int(geomSize.size()); i++)
1154 {
1155 switch(geomSize[i])
1156 {
1157 case 1:
1158 fout << "1 "; //Point
1159 break;
1160 case 2:
1161 fout << "3 "; //Line
1162 break;
1163 case 3:
1164 fout << "5 "; //Triangle
1165 break;
1166 default:
1167 fout << "7 "; //Polygon
1168 }
1169 if(i % 30 == 29)
1170 fout << std::endl << indent;
1171 }
1172 fout << std::endl;
1173 fout << " </DataArray>" << std::endl;
1174 fout << " </Cells>" << std::endl;
1175 }
1176
1177 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1179 fout << " </Piece>" << std::endl;
1180 fout << " </UnstructuredGrid>" << std::endl;
1181 fout << "</VTKFile>" << std::endl;
1182 }
1183
1184
1185 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1186 void VisualizationHelpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::writePVTU(std::ofstream& pvtu, std::string baseFname, int numProcs, bool bFineEdges, bool /* bCoarseEdges */) const {
1187 //If using vtk, filenameToWrite now contains final, correct ***.vtu filename (for the current proc)
1188 //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
1189 //pvtu file.
1190 pvtu << "<VTKFile type=\"PUnstructuredGrid\" byte_order=\"LittleEndian\">" << std::endl;
1191 pvtu << " <PUnstructuredGrid GhostLevel=\"0\">" << std::endl;
1192 pvtu << " <PPointData Scalars=\"Node Aggregate Processor\">" << std::endl;
1193 pvtu << " <PDataArray type=\"Int32\" Name=\"Node\"/>" << std::endl;
1194 pvtu << " <PDataArray type=\"Int32\" Name=\"Aggregate\"/>" << std::endl;
1195 pvtu << " <PDataArray type=\"Int32\" Name=\"Processor\"/>" << std::endl;
1196 pvtu << " </PPointData>" << std::endl;
1197 pvtu << " <PPoints>" << std::endl;
1198 pvtu << " <PDataArray type=\"Float64\" NumberOfComponents=\"3\"/>" << std::endl;
1199 pvtu << " </PPoints>" << std::endl;
1200 for(int i = 0; i < numProcs; i++) {
1201 //specify the piece for each proc (the replaceAll expression matches what the filenames will be on other procs)
1202 pvtu << " <Piece Source=\"" << replaceAll(baseFname, "%PROCID", toString(i)) << "\"/>" << std::endl;
1203 }
1204 //reference files with graph pieces, if applicable
1205 if(bFineEdges)
1206 {
1207 for(int i = 0; i < numProcs; i++)
1208 {
1209 std::string fn = replaceAll(baseFname, "%PROCID", toString(i));
1210 pvtu << " <Piece Source=\"" << fn.insert(fn.rfind(".vtu"), "-finegraph") << "\"/>" << std::endl;
1211 }
1212 }
1213 /*if(doCoarseGraphEdges_)
1214 {
1215 for(int i = 0; i < numProcs; i++)
1216 {
1217 std::string fn = replaceAll(baseFname, "%PROCID", toString(i));
1218 pvtu << " <Piece Source=\"" << fn.insert(fn.rfind(".vtu"), "-coarsegraph") << "\"/>" << std::endl;
1219 }
1220 }*/
1221 pvtu << " </PUnstructuredGrid>" << std::endl;
1222 pvtu << "</VTKFile>" << std::endl;
1223 pvtu.close();
1224 }
1225
1226 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1228 try {
1229 std::ofstream color("random-colormap.xml");
1230 color << "<ColorMap name=\"MueLu-Random\" space=\"RGB\">" << std::endl;
1231 //Give -1, -2, -3 distinctive colors (so that the style functions can have constrasted geometry types)
1232 //Do red, orange, yellow to constrast with cool color spectrum for other types
1233 color << " <Point x=\"" << -1 << "\" o=\"1\" r=\"1\" g=\"0\" b=\"0\"/>" << std::endl;
1234 color << " <Point x=\"" << -2 << "\" o=\"1\" r=\"1\" g=\"0.6\" b=\"0\"/>" << std::endl;
1235 color << " <Point x=\"" << -3 << "\" o=\"1\" r=\"1\" g=\"1\" b=\"0\"/>" << std::endl;
1236 srand(time(NULL));
1237 for(int i = 0; i < 5000; i += 4) {
1238 color << " <Point x=\"" << i << "\" o=\"1\" r=\"" << (rand() % 50) / 256.0 << "\" g=\"" << (rand() % 256) / 256.0 << "\" b=\"" << (rand() % 256) / 256.0 << "\"/>" << std::endl;
1239 }
1240 color << "</ColorMap>" << std::endl;
1241 color.close();
1242 }
1243 catch(std::exception& e) {
1244 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError,
1245 "VisualizationHelpers::buildColormap: Error while building colormap file: " << e.what());
1246 }
1247 }
1248
1249} // namespace MueLu
1250
1251#endif /* MUELU_VISUALIZATIONHELPERS_DEF_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception throws to report errors in the internal logical of the program.
void writeFileVTKData(std::ofstream &fout, std::vector< int > &uniqueFine, LocalOrdinal myAggOffset, ArrayRCP< LocalOrdinal > &vertex2AggId, int myRank) const
static double pointDistFromTri(myVec3 point, myVec3 v1, myVec3 v2, myVec3 v3)
void writeFileVTKClosing(std::ofstream &fout) const
static myVec3 crossProduct(myVec3 v1, myVec3 v2)
std::string replaceAll(std::string result, const std::string &replaceWhat, const std::string &replaceWithWhat) const
static myVec2 vecSubtract(myVec2 v1, myVec2 v2)
std::string getBaseFileName(int numProcs, int level, const Teuchos::ParameterList &pL) const
std::string getFileName(int numProcs, int myRank, int level, const Teuchos::ParameterList &pL) const
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
std::vector< int > makeUnique(std::vector< int > &vertices) const
replaces node indices in vertices with compressed unique indices, and returns list of unique points
static double dotProduct(myVec2 v1, myVec2 v2)
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 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)
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)
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
static bool isInFront(myVec3 point, myVec3 inPlane, myVec3 n)
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 void doPointCloud(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes)
RCP< ParameterList > GetValidParameterList() const
void writePVTU(std::ofstream &pvtu, std::string baseFname, int numProcs, bool bFineEdges=false, bool bCoarseEdges=false) const
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)
static void doJacks(std::vector< int > &vertices, std::vector< int > &geomSizes, LO numLocalAggs, LO numFineNodes, const std::vector< bool > &isRoot, const ArrayRCP< LO > &vertex2AggId)
void writeFileVTKNodes(std::ofstream &fout, std::vector< int > &uniqueFine, Teuchos::RCP< const Map > &nodeMap) const
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)
Namespace for MueLu classes and methods.
void replaceAll(std::string &str, const std::string &from, const std::string &to)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.