LBIBCell
 All Classes Functions Variables Friends Pages
vtkForceReporter.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/vtkForceReporter.hpp>
23 #include <UtilLib/include/Exception.hpp>
24 #include <LbmLib/include/nodes/PhysicalNode.hpp>
25 #include <LbmLib/include/solver/ForceSolver.hpp>
26 #include <LbmLib/include/solver/AbstractForceStruct.hpp>
27 
28 #include <vtkSmartPointer.h>
29 #include <vtkPolyData.h>
30 #include <vtkLine.h>
31 #include <vtkCellData.h>
32 #include <vtkCellArray.h>
33 #include <vtkUnsignedIntArray.h>
34 #include <vtkMultiBlockDataSet.h>
35 #include <vtkXMLMultiBlockDataWriter.h>
36 #include <vtkXMLMultiBlockDataReader.h>
37 
38 #include <sys/stat.h>
39 #include <sstream>
40 #include <fstream>
41 #include <iomanip>
42 
43 namespace LbmLib {
44 namespace reportHandler {
45 void vtkForceReporter::operator()(unsigned int time) const {
46  std::stringstream filename;
47  filename << filename_ << "_" << time << ".vtm";
48 
49  vtkSmartPointer<vtkMultiBlockDataSet> multiBDS =
50  vtkSmartPointer<vtkMultiBlockDataSet>::New (); // multiblock container
51  vtkSmartPointer<vtkXMLMultiBlockDataReader> reader =
52  vtkSmartPointer<vtkXMLMultiBlockDataReader>::New(); // reader
53  vtkSmartPointer<vtkXMLMultiBlockDataWriter> writer =
54  vtkSmartPointer<vtkXMLMultiBlockDataWriter>::New(); // writer
55  vtkSmartPointer<vtkPolyData> linePolyData =
56  vtkSmartPointer<vtkPolyData>::New(); // line container
57  vtkSmartPointer<vtkPoints> forcepoints =
58  vtkSmartPointer<vtkPoints>::New();
59  vtkSmartPointer<vtkCellArray> lines =
60  vtkSmartPointer<vtkCellArray>::New();
61  std::vector<vtkSmartPointer<vtkLine> > linevector;
62  vtkSmartPointer<vtkUnsignedIntArray> forcetype =
63  vtkSmartPointer<vtkUnsignedIntArray>::New();
64  forcetype->SetNumberOfComponents(1);
65  forcetype->SetName("ForceType");
66 
67  LbmLib::solver::map_forcestruct forcemap = this->forcesolver_.getAllForces();
68  unsigned int pointcounter = 0;
69 
70  for (auto iterator=forcemap.begin();
71  iterator != forcemap.end();
72  iterator++) {
73  for (int k=0; k<iterator->second.size();++k) {
74  if (iterator->second[k]->getPartnerGeometryNode() != 0) { // only if there is a line to draw
75  forcepoints->InsertNextPoint(this->geometryNodes_.at(iterator->first)->getXPos(),
76  this->geometryNodes_.at(iterator->first)->getYPos(),
77  0.0);
78  forcepoints->InsertNextPoint(this->geometryNodes_.at(iterator->second[k]->getPartnerGeometryNode())->getXPos(),
79  this->geometryNodes_.at(iterator->second[k]->getPartnerGeometryNode())->getYPos(),
80  0.0);
81  forcetype->InsertNextValue(iterator->second[k]->getType()); // force type
82  linevector.push_back(vtkSmartPointer<vtkLine>::New()); // a new line
83  linevector.back()->GetPointIds()->SetId(0,pointcounter); // first point of the line
84  linevector.back()->GetPointIds()->SetId(1,pointcounter+1); // second point of the line
85  pointcounter += 2; // increment because two points have been added
86  }
87  }
88  }
89 
90  for (int k=0; k<linevector.size(); ++k) {
91  lines->InsertNextCell(linevector[k]);
92  }
93 
94  //add the points and lines to the polydata
95  linePolyData->SetPoints(forcepoints);
96  linePolyData->SetLines(lines);
97  linePolyData->GetCellData()->AddArray(forcetype);
98 
99  // assembly of the multiblock container:
100  reader->SetFileName(filename.str().c_str());
101  struct stat buffer;
102  if (stat(filename.str().c_str(), &buffer) == 0) { // check if file already existing
103  reader->Update(); // if yes, load its content
104  multiBDS->ShallowCopy(reader->GetOutput());
105  }
106 
107  multiBDS->SetBlock(multiBDS->GetNumberOfBlocks(),linePolyData); // add the force lines
108 
109  writer->SetFileName(filename.str().c_str());
110 #if VTK_MAJOR_VERSION <= 5
111  writer->SetInput(multiBDS);
112 #else
113  writer->SetInputData(multiBDS);
114 #endif
115  //writer->SetCompressorTypeToNone();
116  writer->SetDataModeToAscii();
117  //writer->SetDataModeToBinary();
118  writer->Write();
119 }
120 }
121 } // end namespace
virtual void operator()(unsigned int time) const
operator() Writes the report
const map_forcestruct getAllForces() const
getAllForces returns the entire map {nodeid, AbstractForceStruct}
const std::string filename_
filename_ Stores the filename of this functor