LBIBCell
 All Classes Functions Variables Friends Pages
vtkCellReporter.cpp
1 /* Copyright (c) 2013 David Sichau <mail"at"sichau"dot"eu>
2  * 2013-2015 Simon Tanaka <tanakas"at"gmx"dot"ch>
3  *
4  * Permission is hereby granted, free of charge, to any person obtaining a copy
5  * of this software and associated documentation files (the "Software"), to deal
6  * in the Software without restriction, including without limitation the rights
7  * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
8  * copies of the Software, and to permit persons to whom the Software is
9  * furnished to do so, subject to the following conditions:
10  *
11  * The above copyright notice and this permission notice shall be included in
12  * all copies or substantial portions of the Software.
13  *
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND EXPRESS OR
15  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
16  * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
17  * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
18  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
19  * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
20  * THE SOFTWARE.
21  */
22 #include <LbmLib/include/reportHandler/vtkCellReporter.hpp>
23 #include <UtilLib/include/Exception.hpp>
24 
25 #include <vtkSmartPointer.h>
26 #include <vtkPolygon.h>
27 #include <vtkCellArray.h>
28 #include <vtkPolyData.h>
29 #include <vtkDoubleArray.h>
30 #include <vtkStringArray.h>
31 #include <vtkPointData.h>
32 #include <vtkMultiBlockDataSet.h>
33 #include <vtkXMLMultiBlockDataWriter.h>
34 #include <vtkXMLMultiBlockDataReader.h>
35 
36 #include <sys/stat.h>
37 #include <sstream>
38 #include <fstream>
39 #include <iomanip>
40 
41 namespace LbmLib {
42 namespace reportHandler {
43 void vtkCellReporter::operator()(unsigned int time) const {
44  std::stringstream filename;
45  filename << filename_ << "_" << time << ".vtm";
46 
47  vtkSmartPointer<vtkPoints> polygonpoints =
48  vtkSmartPointer<vtkPoints>::New();
49  vtkSmartPointer<vtkCellArray> polygons =
50  vtkSmartPointer<vtkCellArray>::New();
51  vtkSmartPointer<vtkPolyData> polygonPolyData =
52  vtkSmartPointer<vtkPolyData>::New();
53  vtkSmartPointer<vtkMultiBlockDataSet> multiBDS =
54  vtkSmartPointer<vtkMultiBlockDataSet>::New ();
55  vtkSmartPointer<vtkXMLMultiBlockDataReader> reader =
56  vtkSmartPointer<vtkXMLMultiBlockDataReader>::New();
57  vtkSmartPointer<vtkXMLMultiBlockDataWriter> writer =
58  vtkSmartPointer<vtkXMLMultiBlockDataWriter>::New();
59  std::vector<vtkSmartPointer<vtkPolygon> > polygonvector;
60  std::map<unsigned int,std::vector<std::shared_ptr<LbmLib::geometry::Connection> > > celldefinition;
61  std::shared_ptr<LbmLib::geometry::Connection> startC = nullptr;
62  std::shared_ptr<LbmLib::geometry::Connection> tempC = nullptr;
63  vtkSmartPointer<vtkDoubleArray> domainidentifier =
64  vtkSmartPointer<vtkDoubleArray>::New();
65  domainidentifier->SetName("domain_identifier");
66  vtkSmartPointer<vtkDoubleArray> celltype =
67  vtkSmartPointer<vtkDoubleArray>::New();
68  celltype->SetName("cell_type");
69 
70  vtkSmartPointer<vtkStringArray> stringCDESolverArray =
71  vtkSmartPointer<vtkStringArray>::New();
72 
73  stringCDESolverArray->SetNumberOfComponents(1);
74  stringCDESolverArray->SetName("boundarydescriptors");
75 
76  std::stringstream boundarydescriptor;
77 
78  int globalpointcounter = 0;
79  int polygonpointcounter;
80  int polygoncounter = 0;
81  unsigned int temp_celltype;
82 
83  for (auto it : this->connections_) { // create a cell definition map:
84  celldefinition[(*it).getDomainIdentifier()].push_back(it);
85  }
86 
87  for (auto it : celldefinition) { // loop the cells
88  startC = it.second[0];
89  tempC = startC;
90  polygonpointcounter = 0;
91  polygonvector.push_back(vtkSmartPointer<vtkPolygon>::New());
92 
93  do {
94  polygonpoints->InsertNextPoint(tempC->getGeometryNodes().second->getXPos(),
95  tempC->getGeometryNodes().second->getYPos(),
96  0.0); // add the points
97  domainidentifier->InsertNextValue(tempC->getDomainIdentifier()); // add the domain identifier
98 
99  if ( this->cellTypeTrackerMap_.find(tempC->getDomainIdentifier()) == this->cellTypeTrackerMap_.end() ) {
100  lbm_fail("domainIdentifier not found in the cellTrackerMap.");
101  } else {
102  temp_celltype = this->cellTypeTrackerMap_.at(tempC->getDomainIdentifier()); // get the celltype
103  }
104 
105  celltype->InsertNextValue(temp_celltype); // add the cell type
106 
107  boundarydescriptor.str(""); // clear string
108  boundarydescriptor << "\t"; // to avoid segfaults when there is no other content
109  for (const auto& bSolver : tempC->getBoundaryConditionDescriptor()) { // assmble a string with the boundary descriptor
110  for (const auto& cdeSolver : bSolver.second) {
111  boundarydescriptor << bSolver.first << '\t' << cdeSolver << '\t';
112  }
113  }
114 
115  stringCDESolverArray->InsertNextValue(boundarydescriptor.str().c_str());
116  polygonvector.back()->GetPointIds()->InsertId(polygonpointcounter, globalpointcounter);
117 
118  tempC = tempC->getGeometryNodes().second->getConnection<1>(); // go to next connection
119  globalpointcounter++;
120  polygonpointcounter++;
121  } while(tempC != startC); // do until reaching beginning of the polygon
122  polygoncounter++;
123  }
124 
125  for (auto k : polygonvector) { // add the polygons to a list of polygons
126  polygons->InsertNextCell(k);
127  }
128 
129  polygonPolyData->SetPoints(polygonpoints); // add the points to the polydata container
130  polygonPolyData->SetPolys(polygons);
131  polygonPolyData->GetPointData()->AddArray(domainidentifier); // add the attributes
132  polygonPolyData->GetPointData()->AddArray(celltype); // add the attributes
133  polygonPolyData->GetPointData()->AddArray(stringCDESolverArray); // add the boundarydescriptor strings
134 
135  reader->SetFileName(filename.str().c_str());
136  struct stat buffer;
137  if (stat(filename.str().c_str(), &buffer) == 0) { // check if file already existing
138 
139  reader->Update(); // if yes, load its content
140  multiBDS->DeepCopy(reader->GetOutput());
141  }
142 
143  multiBDS->SetBlock(multiBDS->GetNumberOfBlocks(),polygonPolyData);
144 
145  writer->SetFileName(filename.str().c_str());
146 #if VTK_MAJOR_VERSION <= 5
147  writer->SetInput(multiBDS);
148 #else
149  writer->SetInputData(multiBDS);
150 #endif
151  //writer->SetCompressorTypeToNone();
152  writer->SetDataModeToAscii();
153  //writer->SetDataModeToBinary();
154  writer->Write();
155 }
156 }
157 } // end namespace
virtual void operator()(unsigned int time) const
operator() Writes the report
const std::string filename_
filename_ Stores the filename of this functor