53 #ifndef MUELU_AGGREGATIONEXPORTFACTORY_DEF_HPP_
54 #define MUELU_AGGREGATIONEXPORTFACTORY_DEF_HPP_
56 #include <Xpetra_Matrix.hpp>
57 #include <Xpetra_CrsMatrixWrap.hpp>
58 #include <Xpetra_ImportFactory.hpp>
59 #include <Xpetra_MultiVectorFactory.hpp>
62 #include "MueLu_Aggregates.hpp"
63 #include "MueLu_Graph.hpp"
64 #include "MueLu_AmalgamationFactory.hpp"
65 #include "MueLu_AmalgamationInfo.hpp"
67 #include "MueLu_Utilities.hpp"
76 #ifdef HAVE_MUELU_CGAL //Include all headers needed for both 2D and 3D fixed-alpha alpha shapes
77 #include "CGAL/Exact_predicates_inexact_constructions_kernel.h"
78 #include "CGAL/Delaunay_triangulation_2.h"
79 #include "CGAL/Delaunay_triangulation_3.h"
80 #include "CGAL/Alpha_shape_2.h"
81 #include "CGAL/Fixed_alpha_shape_3.h"
82 #include "CGAL/algorithm.h"
87 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
91 std::string output_msg =
"Output filename template (%TIMESTEP is replaced by \'Output file: time step\' variable,"
92 "%ITER is replaced by \'Output file: iter\' variable, %LEVELID is replaced level id, %PROCID is replaced by processor id)";
93 std::string output_def =
"aggs_level%LEVELID_proc%PROCID.out";
102 validParamList->
set< std::string > (
"Output filename", output_def, output_msg);
103 validParamList->
set<
int > (
"Output file: time step", 0,
"time step variable for output file name");
104 validParamList->
set<
int > (
"Output file: iter", 0,
"nonlinear iteration variable for output file name");
107 validParamList->
set< std::string > (
"aggregation: output filename",
"",
"filename for VTK-style visualization output");
108 validParamList->
set<
int > (
"aggregation: output file: time step", 0,
"time step variable for output file name");
109 validParamList->
set<
int > (
"aggregation: output file: iter", 0,
"nonlinear iteration variable for output file name");
110 validParamList->
set<std::string> (
"aggregation: output file: agg style",
"Point Cloud",
"style of aggregate visualization for VTK output");
111 validParamList->
set<
bool> (
"aggregation: output file: fine graph edges",
false,
"Whether to draw all fine node connections along with the aggregates.");
112 validParamList->
set<
bool> (
"aggregation: output file: coarse graph edges",
false,
"Whether to draw all coarse node connections along with the aggregates.");
113 validParamList->
set<
bool> (
"aggregation: output file: build colormap",
false,
"Whether to output a random colormap for ParaView in a separate XML file.");
114 return validParamList;
117 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
119 Input(fineLevel,
"Aggregates");
120 Input(fineLevel,
"DofsPerNode");
121 Input(fineLevel,
"UnAmalgamationInfo");
125 if(pL.
isParameter(
"aggregation: output filename") && pL.
get<std::string>(
"aggregation: output filename").length())
127 Input(fineLevel,
"Coordinates");
128 Input(fineLevel,
"A");
129 Input(fineLevel,
"Graph");
130 if(pL.
get<
bool>(
"aggregation: output file: coarse graph edges"))
132 Input(coarseLevel,
"Coordinates");
133 Input(coarseLevel,
"A");
134 Input(coarseLevel,
"Graph");
139 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147 int numProcs = comm->getSize();
148 int myRank = comm->getRank();
149 string masterFilename = pL.
get<std::string>(
"aggregation: output filename");
150 string pvtuFilename =
"";
151 string localFilename = pL.get<std::string>(
"Output filename");
152 string filenameToWrite;
154 doCoarseGraphEdges_ = pL.get<
bool>(
"aggregation: output file: coarse graph edges");
155 doFineGraphEdges_ = pL.get<
bool>(
"aggregation: output file: fine graph edges");
156 if(masterFilename.length())
159 filenameToWrite = masterFilename;
160 if(filenameToWrite.rfind(
".vtu") == string::npos)
161 filenameToWrite.append(
".vtu");
162 if(numProcs > 1 && filenameToWrite.rfind(
"%PROCID") == string::npos)
163 filenameToWrite.insert(filenameToWrite.rfind(
".vtu"),
"-proc%PROCID");
166 filenameToWrite = localFilename;
167 LocalOrdinal DofsPerNode = Get< LocalOrdinal > (fineLevel,
"DofsPerNode");
171 if(doCoarseGraphEdges_)
172 Ac = Get<RCP<Matrix> >(coarseLevel,
"A");
177 if(doFineGraphEdges_)
178 fineGraph = Get<RCP<GraphBase> >(fineLevel,
"Graph");
179 if(doCoarseGraphEdges_)
180 coarseGraph = Get<RCP<GraphBase> >(coarseLevel,
"Graph");
183 coords = Get<RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,
LocalOrdinal,
GlobalOrdinal,
Node> > >(fineLevel,
"Coordinates");
184 if(doCoarseGraphEdges_)
185 coordsCoarse = Get<RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,
LocalOrdinal,
GlobalOrdinal,
Node> > >(coarseLevel,
"Coordinates");
186 dims_ = coords->getNumVectors();
190 RCP<Import> coordImporter = Xpetra::ImportFactory<LocalOrdinal, GlobalOrdinal, Node>::Build(coords->getMap(), Amat->getColMap());
192 ghostedCoords->doImport(*coords, *coordImporter, Xpetra::INSERT);
193 coords = ghostedCoords;
195 if(doCoarseGraphEdges_)
197 RCP<Import> coordImporter = Xpetra::ImportFactory<LocalOrdinal, GlobalOrdinal, Node>::Build(coordsCoarse->getMap(), Ac->getColMap());
199 ghostedCoords->doImport(*coordsCoarse, *coordImporter, Xpetra::INSERT);
200 coordsCoarse = ghostedCoords;
204 GetOStream(
Runtime0) <<
"AggregationExportFactory: DofsPerNode: " << DofsPerNode << std::endl;
210 vertex2AggId_ = vertex2AggId;
213 std::vector<GlobalOrdinal> numAggsGlobal (numProcs, 0);
214 std::vector<GlobalOrdinal> numAggsLocal (numProcs, 0);
215 std::vector<GlobalOrdinal> minGlobalAggId(numProcs, 0);
217 numAggsLocal[myRank] = aggregates->GetNumAggregates();
218 Teuchos::reduceAll(*comm,
Teuchos::REDUCE_SUM, numProcs, &numAggsLocal[0], &numAggsGlobal[0]);
219 for (
int i = 1; i < Teuchos::as<int>(numAggsGlobal.size()); ++i)
221 numAggsGlobal [i] += numAggsGlobal[i-1];
222 minGlobalAggId[i] = numAggsGlobal[i-1];
227 aggsOffset_ = minGlobalAggId[myRank];
230 amalgInfo->UnamalgamateAggregates(*aggregates, aggStart, aggToRowMap);
231 int timeStep = pL.
get<
int > (
"Output file: time step");
232 int iter = pL.get<
int > (
"Output file: iter");
238 string masterStem =
"";
241 masterStem = filenameToWrite.substr(0, filenameToWrite.rfind(
".vtu"));
242 masterStem = this->
replaceAll(masterStem,
"%PROCID",
"");
244 pvtuFilename = masterStem +
"-master.pvtu";
245 string baseFname = filenameToWrite;
247 GetOStream(
Runtime0) <<
"AggregationExportFactory: outputfile \"" << filenameToWrite <<
"\"" << std::endl;
248 ofstream fout(filenameToWrite.c_str());
249 GO numAggs = aggregates->GetNumAggregates();
252 GO indexBase = aggregates->GetMap()->getIndexBase();
253 GO offset = amalgInfo->GlobalOffset();
254 vector<GlobalOrdinal> nodeIds;
255 for (
int i = 0; i < numAggs; ++i) {
256 fout <<
"Agg " << minGlobalAggId[myRank] + i <<
" Proc " << myRank <<
":";
259 for (
int k = aggStart[i]; k < aggStart[i+1]; ++k) {
260 nodeIds.push_back((aggToRowMap[k] - offset - indexBase) / DofsPerNode + indexBase);
264 std::sort(nodeIds.begin(), nodeIds.end());
265 typename std::vector<GlobalOrdinal>::iterator endLocation = std::unique(nodeIds.begin(), nodeIds.end());
266 nodeIds.erase(endLocation, nodeIds.end());
269 for(
typename std::vector<GlobalOrdinal>::iterator printIt = nodeIds.begin(); printIt != nodeIds.end(); printIt++)
270 fout <<
" " << *printIt;
281 numNodes_ = coords->getLocalLength();
285 zCoords_ = Teuchos::null;
286 if(doCoarseGraphEdges_)
295 if(doCoarseGraphEdges_)
299 aggSizes_ = aggregates->ComputeAggregateSizes();
301 string aggStyle =
"Point Cloud";
304 aggStyle = pL.get<
string>(
"aggregation: output file: agg style");
306 catch(std::exception& e) {}
307 vector<int> vertices;
308 vector<int> geomSizes;
310 nodeMap_ = Amat->getMap();
313 isRoot_.push_back(aggregates->IsRoot(i));
317 if(myRank == 0 && (numProcs != 1 || doCoarseGraphEdges_ || doFineGraphEdges_))
319 ofstream pvtu(pvtuFilename.c_str());
320 writePVTU_(pvtu, baseFname, numProcs);
323 if(aggStyle ==
"Point Cloud")
324 this->doPointCloud(vertices, geomSizes, numAggs_, numNodes_);
325 else if(aggStyle ==
"Jacks")
326 this->doJacks(vertices, geomSizes, numAggs_, numNodes_, isRoot_, vertex2AggId_);
327 else if(aggStyle ==
"Jacks++")
328 doJacksPlus_(vertices, geomSizes);
329 else if(aggStyle ==
"Convex Hulls")
330 doConvexHulls(vertices, geomSizes);
331 else if(aggStyle ==
"Alpha Hulls")
333 #ifdef HAVE_MUELU_CGAL
334 doAlphaHulls_(vertices, geomSizes);
336 GetOStream(
Warnings0) <<
" Trilinos was not configured with CGAL so Alpha Hulls not available.\n Using Convex Hulls instead." << std::endl;
337 doConvexHulls(vertices, geomSizes);
342 GetOStream(
Warnings0) <<
" Unrecognized agg style.\n Possible values are Point Cloud, Jacks, Jacks++, Convex Hulls and Alpha Hulls.\n Defaulting to Point Cloud." << std::endl;
343 aggStyle =
"Point Cloud";
344 this->doPointCloud(vertices, geomSizes, numAggs_, numNodes_);
346 writeFile_(fout, aggStyle, vertices, geomSizes);
347 if(doCoarseGraphEdges_)
349 string fname = filenameToWrite;
350 string cEdgeFile = fname.insert(fname.rfind(
".vtu"),
"-coarsegraph");
351 std::ofstream edgeStream(cEdgeFile.c_str());
352 doGraphEdges_(edgeStream, Ac, coarseGraph,
false, DofsPerNode);
354 if(doFineGraphEdges_)
356 string fname = filenameToWrite;
357 string fEdgeFile = fname.insert(fname.rfind(
".vtu"),
"-finegraph");
358 std::ofstream edgeStream(fEdgeFile.c_str());
359 doGraphEdges_(edgeStream, Amat, fineGraph,
true, DofsPerNode);
361 if(myRank == 0 && pL.get<
bool>(
"aggregation: output file: build colormap"))
369 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
375 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
379 this->doConvexHulls2D(vertices, geomSizes, numAggs_, numNodes_, isRoot_, vertex2AggId_, xCoords_, yCoords_);
381 this->doConvexHulls3D(vertices, geomSizes, numAggs_, numNodes_, isRoot_, vertex2AggId_, xCoords_, yCoords_, zCoords_);
384 #ifdef HAVE_MUELU_CGAL
385 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
390 doAlphaHulls2D_(vertices, geomSizes);
392 doAlphaHulls3D_(vertices, geomSizes);
395 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
401 typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
403 typedef K::Point_2 Point;
404 typedef K::Segment_2 Segment;
405 typedef CGAL::Alpha_shape_vertex_base_2<K> Vb;
406 typedef CGAL::Alpha_shape_face_base_2<K> Fb;
407 typedef CGAL::Triangulation_data_structure_2<Vb,Fb> Tds;
408 typedef CGAL::Delaunay_triangulation_2<K,Tds> Triangulation_2;
409 typedef CGAL::Alpha_shape_2<Triangulation_2> Alpha_shape_2;
410 typedef Alpha_shape_2::Alpha_shape_edges_iterator Alpha_shape_edges_iterator;
411 #if 0 // taw: does not compile with CGAL 4.8
412 for(
int i = 0; i < numAggs_; i++)
415 list<Point> aggPoints;
416 vector<int> aggNodes;
417 for(
int j = 0; j < numNodes_; j++)
419 if(vertex2AggId_[j] == i)
421 Point p(xCoords_[j], yCoords_[j]);
422 aggPoints.push_back(p);
423 aggNodes.push_back(j);
426 Alpha_shape_2 hull(aggPoints.begin(), aggPoints.end(), FT(ALPHA_VAL), Alpha_shape_2::GENERAL);
427 vector<Segment> segments;
428 CGAL::alpha_edges(hull, back_inserter(segments));
429 vertices.reserve(vertices.size() + 2 * segments.size());
430 geomSizes.reserve(geomSizes.size() + segments.size());
431 for(
size_t j = 0; j < segments.size(); j++)
433 for(
size_t k = 0; k < aggNodes.size(); k++)
435 if(fabs(segments[j][0].x == xCoords_[aggNodes[k]]) < 1e-12 && fabs(segments[j][0].y == yCoords_[aggNodes[k]]) < 1e-12)
437 vertices.push_back(aggNodes[k]);
441 for(
size_t k = 0; k < aggNodes.size(); k++)
443 if(fabs(segments[j][1].x == xCoords_[aggNodes[k]]) < 1e-12 && fabs(segments[j][1].y == yCoords_[aggNodes[k]]) < 1e-12)
445 vertices.push_back(aggNodes[k]);
449 geomSizes.push_back(2);
455 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
458 typedef CGAL::Exact_predicates_inexact_constructions_kernel Gt;
459 #if 0 // does not compile with CGAL 4-8
460 typedef CGAL::Alpha_shape_cell_base_3<Gt> Fb;
461 typedef CGAL::Triangulation_data_structure_3<Vb,Fb> Tds;
462 typedef CGAL::Delaunay_triangulation_3<Gt,Tds> Triangulation_3;
463 typedef Gt::Point_3 Point;
464 typedef Alpha_shape_3::Alpha_iterator Alpha_iterator;
465 typedef Alpha_shape_3::Cell_handle Cell_handle;
466 typedef Alpha_shape_3::Vertex_handle Vertex_handle;
467 typedef Alpha_shape_3::Facet Facet;
468 typedef Alpha_shape_3::Edge Edge;
469 typedef Gt::Weighted_point Weighted_point;
470 typedef Gt::Bare_point Bare_point;
471 const double ALPHA_VAL = 2;
474 for(
int i = 0; i < numAggs_; i++)
476 list<Point> aggPoints;
477 vector<int> aggNodes;
478 for(
int j = 0; j < numNodes_; j++)
480 if(vertex2AggId[j] == i)
482 Point p(xCoords_[j], yCoords_[j], zCoords_[j]);
483 aggPoints.push_back(p);
484 aggNodes.push_back(j);
487 Fixed_alpha_shape_3 hull(aggPoints.begin(), aggPoints.end(), FT(ALPHA_VAL));
488 list<Cell_handle> cells;
491 hull.get_alpha_shape_cells(back_inserter(cells));
492 hull.get_alpha_shape_facets(back_inserter(facets));
493 hull.get_alpha_shape_edges(back_inserter(edges));
494 for(
size_t j = 0; j < cells.size(); j++)
497 tetPoints[0] = cells[j]->vertex(0);
498 tetPoints[1] = cells[j]->vertex(1);
499 tetPoints[2] = cells[j]->vertex(2);
500 tetPoints[3] = cells[j]->vertex(3);
501 for(
int k = 0; k < 4; k++)
503 for(
size_t l = 0; l < aggNodes.size(); l++)
505 if(fabs(tetPoints[k].x - xCoords_[aggNodes[l]]) < 1e-12 &&
506 fabs(tetPoints[k].y - yCoords_[aggNodes[l]]) < 1e-12 &&
507 fabs(tetPoints[k].z - zCoords_[aggNodes[l]]) < 1e-12)
509 vertices.push_back(aggNodes[l]);
514 geomSizes.push_back(-10);
516 for(
size_t j = 0; j < facets.size(); j++)
519 indices[0] = (facets[i].second + 1) % 4;
520 indices[1] = (facets[i].second + 2) % 4;
521 indices[2] = (facets[i].second + 3) % 4;
522 if(facets[i].second % 2 == 0)
523 swap(indices[0], indices[1]);
525 facetPts[0] = facets[i].first->vertex(indices[0])->point();
526 facetPts[1] = facets[i].first->vertex(indices[1])->point();
527 facetPts[2] = facets[i].first->vertex(indices[2])->point();
529 for(
size_t l = 0; l < aggNodes.size(); l++)
531 if(fabs(facetPts[k].x - xCoords_[aggNodes[l]]) < 1e-12 &&
532 fabs(facetPts[k].y - yCoords_[aggNodes[l]]) < 1e-12 &&
533 fabs(facetPts[k].z - zCoords_[aggNodes[l]]) < 1e-12)
535 vertices.push_back(aggNodes[l]);
539 geomSizes.push_back(3);
541 for(
size_t j = 0; j < edges.size(); j++)
550 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
557 vector<pair<int, int> > vert1;
558 vector<pair<int, int> > vert2;
559 if(A->isGloballyIndexed())
565 A->getGlobalRowView(globRow, indices, values);
566 neighbors = G->getNeighborVertices((
LocalOrdinal) globRow);
569 while(gEdge !=
int(neighbors.size()))
573 if(neighbors[gEdge] == indices[aEdge])
576 vert1.push_back(pair<int, int>(
int(globRow), neighbors[gEdge]));
583 vert2.push_back(pair<int, int>(
int(globRow), neighbors[gEdge]));
589 vert1.push_back(pair<int, int>(
int(globRow), neighbors[gEdge]));
601 A->getLocalRowView(locRow, indices, values);
602 neighbors = G->getNeighborVertices(locRow);
606 while(gEdge !=
int(neighbors.size()))
610 if(neighbors[gEdge] == indices[aEdge])
612 vert1.push_back(pair<int, int>(locRow, neighbors[gEdge]));
618 vert2.push_back(pair<int, int>(locRow, neighbors[gEdge]));
624 vert1.push_back(pair<int, int>(locRow, neighbors[gEdge]));
630 for(
size_t i = 0; i < vert1.size(); i ++)
632 if(vert1[i].first > vert1[i].second)
634 int temp = vert1[i].first;
635 vert1[i].first = vert1[i].second;
636 vert1[i].second = temp;
639 for(
size_t i = 0; i < vert2.size(); i++)
641 if(vert2[i].first > vert2[i].second)
643 int temp = vert2[i].first;
644 vert2[i].first = vert2[i].second;
645 vert2[i].second = temp;
648 sort(vert1.begin(), vert1.end());
649 vector<pair<int, int> >::iterator newEnd = unique(vert1.begin(), vert1.end());
650 vert1.erase(newEnd, vert1.end());
651 sort(vert2.begin(), vert2.end());
652 newEnd = unique(vert2.begin(), vert2.end());
653 vert2.erase(newEnd, vert2.end());
655 points1.reserve(2 * vert1.size());
656 for(
size_t i = 0; i < vert1.size(); i++)
658 points1.push_back(vert1[i].first);
659 points1.push_back(vert1[i].second);
662 points2.reserve(2 * vert2.size());
663 for(
size_t i = 0; i < vert2.size(); i++)
665 points2.push_back(vert2[i].first);
666 points2.push_back(vert2[i].second);
668 vector<int> unique1 = this->makeUnique(points1);
669 vector<int> unique2 = this->makeUnique(points2);
670 fout <<
"<VTKFile type=\"UnstructuredGrid\" byte_order=\"LittleEndian\">" << endl;
671 fout <<
" <UnstructuredGrid>" << endl;
672 fout <<
" <Piece NumberOfPoints=\"" << unique1.size() + unique2.size() <<
"\" NumberOfCells=\"" << vert1.size() + vert2.size() <<
"\">" << endl;
673 fout <<
" <PointData Scalars=\"Node Aggregate Processor\">" << endl;
674 fout <<
" <DataArray type=\"Int32\" Name=\"Node\" format=\"ascii\">" << endl;
677 for(
size_t i = 0; i < unique1.size(); i++)
679 fout << CONTRAST_1_ <<
" ";
681 fout << endl << indent;
683 for(
size_t i = 0; i < unique2.size(); i++)
685 fout << CONTRAST_2_ <<
" ";
686 if((i + 2 * vert1.size()) % 25 == 24)
687 fout << endl << indent;
690 fout <<
" </DataArray>" << endl;
691 fout <<
" <DataArray type=\"Int32\" Name=\"Aggregate\" format=\"ascii\">" << endl;
693 for(
size_t i = 0; i < unique1.size(); i++)
695 fout << CONTRAST_1_ <<
" ";
697 fout << endl << indent;
699 for(
size_t i = 0; i < unique2.size(); i++)
701 fout << CONTRAST_2_ <<
" ";
702 if((i + 2 * vert2.size()) % 25 == 24)
703 fout << endl << indent;
706 fout <<
" </DataArray>" << endl;
707 fout <<
" <DataArray type=\"Int32\" Name=\"Processor\" format=\"ascii\">" << endl;
709 for(
size_t i = 0; i < unique1.size() + unique2.size(); i++)
711 fout << myRank_ <<
" ";
713 fout << endl << indent;
716 fout <<
" </DataArray>" << endl;
717 fout <<
" </PointData>" << endl;
718 fout <<
" <Points>" << endl;
719 fout <<
" <DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">" << endl;
721 for(
size_t i = 0; i < unique1.size(); i++)
725 fout << xCoords_[unique1[i]] <<
" " << yCoords_[unique1[i]] <<
" ";
727 fout << zCoords_[unique1[i]] <<
" ";
731 fout << endl << indent;
735 fout << cx_[unique1[i]] <<
" " << cy_[unique1[i]] <<
" ";
737 fout << cz_[unique1[i]] <<
" ";
741 fout << endl << indent;
744 for(
size_t i = 0; i < unique2.size(); i++)
748 fout << xCoords_[unique2[i]] <<
" " << yCoords_[unique2[i]] <<
" ";
750 fout << zCoords_[unique2[i]] <<
" ";
754 fout << endl << indent;
758 fout << cx_[unique2[i]] <<
" " << cy_[unique2[i]] <<
" ";
760 fout << cz_[unique2[i]] <<
" ";
763 if((i + unique1.size()) % 2)
764 fout << endl << indent;
768 fout <<
" </DataArray>" << endl;
769 fout <<
" </Points>" << endl;
770 fout <<
" <Cells>" << endl;
771 fout <<
" <DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">" << endl;
773 for(
size_t i = 0; i < points1.size(); i++)
775 fout << points1[i] <<
" ";
777 fout << endl << indent;
779 for(
size_t i = 0; i < points2.size(); i++)
781 fout << points2[i] + unique1.size() <<
" ";
782 if((i + points1.size()) % 10 == 9)
783 fout << endl << indent;
786 fout <<
" </DataArray>" << endl;
787 fout <<
" <DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">" << endl;
790 for(
size_t i = 0; i < vert1.size() + vert2.size(); i++)
793 fout << offset <<
" ";
795 fout << endl << indent;
798 fout <<
" </DataArray>" << endl;
799 fout <<
" <DataArray type=\"Int32\" Name=\"types\" format=\"ascii\">" << endl;
801 for(
size_t i = 0; i < vert1.size() + vert2.size(); i++)
805 fout << endl << indent;
808 fout <<
" </DataArray>" << endl;
809 fout <<
" </Cells>" << endl;
810 fout <<
" </Piece>" << endl;
811 fout <<
" </UnstructuredGrid>" << endl;
812 fout <<
"</VTKFile>" << endl;
815 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
819 vector<int> uniqueFine = this->makeUnique(vertices);
821 fout <<
"<!--" << styleName <<
" Aggregates Visualization-->" << endl;
822 fout <<
"<VTKFile type=\"UnstructuredGrid\" byte_order=\"LittleEndian\">" << endl;
823 fout <<
" <UnstructuredGrid>" << endl;
824 fout <<
" <Piece NumberOfPoints=\"" << uniqueFine.size() <<
"\" NumberOfCells=\"" << geomSizes.size() <<
"\">" << endl;
825 fout <<
" <PointData Scalars=\"Node Aggregate Processor\">" << endl;
826 fout <<
" <DataArray type=\"Int32\" Name=\"Node\" format=\"ascii\">" << endl;
830 for(
size_t i = 0; i < uniqueFine.size(); i++)
834 fout << uniqueFine[i] <<
" ";
837 fout << nodeMap_->getGlobalElement(uniqueFine[i]) <<
" ";
839 fout << endl << indent;
842 fout <<
" </DataArray>" << endl;
843 fout <<
" <DataArray type=\"Int32\" Name=\"Aggregate\" format=\"ascii\">" << endl;
845 for(
size_t i = 0; i < uniqueFine.size(); i++)
847 fout << aggsOffset_ + vertex2AggId_[uniqueFine[i]] <<
" ";
849 fout << endl << indent;
852 fout <<
" </DataArray>" << endl;
853 fout <<
" <DataArray type=\"Int32\" Name=\"Processor\" format=\"ascii\">" << endl;
855 for(
size_t i = 0; i < uniqueFine.size(); i++)
857 fout << myRank_ <<
" ";
859 fout << endl << indent;
862 fout <<
" </DataArray>" << endl;
863 fout <<
" </PointData>" << endl;
864 fout <<
" <Points>" << endl;
865 fout <<
" <DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">" << endl;
867 for(
size_t i = 0; i < uniqueFine.size(); i++)
869 fout << xCoords_[uniqueFine[i]] <<
" " << yCoords_[uniqueFine[i]] <<
" ";
873 fout << zCoords_[uniqueFine[i]] <<
" ";
875 fout << endl << indent;
878 fout <<
" </DataArray>" << endl;
879 fout <<
" </Points>" << endl;
880 fout <<
" <Cells>" << endl;
881 fout <<
" <DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">" << endl;
883 for(
size_t i = 0; i < vertices.size(); i++)
885 fout << vertices[i] <<
" ";
887 fout << endl << indent;
890 fout <<
" </DataArray>" << endl;
891 fout <<
" <DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">" << endl;
894 for(
size_t i = 0; i < geomSizes.size(); i++)
896 accum += geomSizes[i];
897 fout << accum <<
" ";
899 fout << endl << indent;
902 fout <<
" </DataArray>" << endl;
903 fout <<
" <DataArray type=\"Int32\" Name=\"types\" format=\"ascii\">" << endl;
905 for(
size_t i = 0; i < geomSizes.size(); i++)
922 fout << endl << indent;
925 fout <<
" </DataArray>" << endl;
926 fout <<
" </Cells>" << endl;
927 fout <<
" </Piece>" << endl;
928 fout <<
" </UnstructuredGrid>" << endl;
929 fout <<
"</VTKFile>" << endl;
932 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
938 ofstream color(
"random-colormap.xml");
939 color <<
"<ColorMap name=\"MueLu-Random\" space=\"RGB\">" << endl;
942 color <<
" <Point x=\"" << CONTRAST_1_ <<
"\" o=\"1\" r=\"1\" g=\"0\" b=\"0\"/>" << endl;
943 color <<
" <Point x=\"" << CONTRAST_2_ <<
"\" o=\"1\" r=\"1\" g=\"0.6\" b=\"0\"/>" << endl;
944 color <<
" <Point x=\"" << CONTRAST_3_ <<
"\" o=\"1\" r=\"1\" g=\"1\" b=\"0\"/>" << endl;
946 for(
int i = 0; i < 5000; i += 4)
948 color <<
" <Point x=\"" << i <<
"\" o=\"1\" r=\"" << (rand() % 50) / 256.0 <<
"\" g=\"" << (rand() % 256) / 256.0 <<
"\" b=\"" << (rand() % 256) / 256.0 <<
"\"/>" << endl;
950 color <<
"</ColorMap>" << endl;
953 catch(std::exception& e)
955 GetOStream(
Warnings0) <<
" Error while building colormap file: " << e.what() << endl;
959 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
966 pvtu <<
"<VTKFile type=\"PUnstructuredGrid\" byte_order=\"LittleEndian\">" << endl;
967 pvtu <<
" <PUnstructuredGrid GhostLevel=\"0\">" << endl;
968 pvtu <<
" <PPointData Scalars=\"Node Aggregate Processor\">" << endl;
969 pvtu <<
" <PDataArray type=\"Int32\" Name=\"Node\"/>" << endl;
970 pvtu <<
" <PDataArray type=\"Int32\" Name=\"Aggregate\"/>" << endl;
971 pvtu <<
" <PDataArray type=\"Int32\" Name=\"Processor\"/>" << endl;
972 pvtu <<
" </PPointData>" << endl;
973 pvtu <<
" <PPoints>" << endl;
974 pvtu <<
" <PDataArray type=\"Float64\" NumberOfComponents=\"3\"/>" << endl;
975 pvtu <<
" </PPoints>" << endl;
976 for(
int i = 0; i < numProcs; i++)
979 pvtu <<
" <Piece Source=\"" << this->
replaceAll(baseFname,
"%PROCID",
toString(i)) <<
"\"/>" << endl;
982 if(doFineGraphEdges_)
984 for(
int i = 0; i < numProcs; i++)
987 pvtu <<
" <Piece Source=\"" << fn.insert(fn.rfind(
".vtu"),
"-finegraph") <<
"\"/>" << endl;
990 if(doCoarseGraphEdges_)
992 for(
int i = 0; i < numProcs; i++)
995 pvtu <<
" <Piece Source=\"" << fn.insert(fn.rfind(
".vtu"),
"-coarsegraph") <<
"\"/>" << endl;
998 pvtu <<
" </PUnstructuredGrid>" << endl;
999 pvtu <<
"</VTKFile>" << endl;
Important warning messages (one line)
void doAlphaHulls3D_(std::vector< int > &vertices, std::vector< int > &geomSizes) const
MueLu::DefaultLocalOrdinal LocalOrdinal
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void doAlphaHulls_(std::vector< int > &vertices, std::vector< int > &geomSizes) const
One-liner description of what is happening.
void doConvexHulls(std::vector< int > &vertices, std::vector< int > &geomSizes) const
void doGraphEdges_(std::ofstream &fout, Teuchos::RCP< Matrix > &A, Teuchos::RCP< GraphBase > &G, bool fine, int dofs) const
void writeFile_(std::ofstream &fout, std::string styleName, std::vector< int > &vertices, std::vector< int > &geomSizes) const
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void writePVTU_(std::ofstream &pvtu, std::string baseFname, int numProcs) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
int GetLevelID() const
Return level number.
Exception throws to report errors in the internal logical of the program.
void doJacksPlus_(std::vector< int > &vertices, std::vector< int > &geomSizes) const
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void replaceAll(std::string &str, const std::string &from, const std::string &to)
void doAlphaHulls2D_(std::vector< int > &vertices, std::vector< int > &geomSizes) const
void buildColormap_() const