LBIBCell
 All Classes Functions Variables Friends Pages
vtkFluidReporter.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/vtkFluidReporter.hpp>
23 #include <UtilLib/include/Exception.hpp>
24 #include <LbmLib/include/nodes/PhysicalNode.hpp>
25 #include <LbmLib/include/solver/CDESolver/CDEAbstractSolver.hpp>
26 #include <LbmLib/include/solver/FluidSolver/FluidSolver.hpp>
27 #include <LbmLib/include/solver/ForceSolver.hpp>
28 
29 #include <vtkSmartPointer.h>
30 #include <vtkDoubleArray.h>
31 #include <vtkPointData.h>
32 #include <vtkMultiBlockDataSet.h>
33 #include <vtkXMLMultiBlockDataWriter.h>
34 #include <vtkXMLMultiBlockDataReader.h>
35 #include <vtkStructuredPoints.h>
36 
37 #include <sys/stat.h>
38 #include <sstream>
39 #include <fstream>
40 #include <iomanip>
41 
42 namespace LbmLib {
43 namespace reportHandler {
44 void vtkFluidReporter::operator()(unsigned int time) const {
45  std::stringstream filename;
46  filename << filename_ << "_" << time << ".vtm";
47 
48  vtkSmartPointer<vtkDoubleArray> densityfield =
49  vtkSmartPointer<vtkDoubleArray>::New();
50  vtkSmartPointer<vtkDoubleArray> velocityfield =
51  vtkSmartPointer<vtkDoubleArray>::New();
52  vtkSmartPointer<vtkStructuredPoints> uniformGrid =
53  vtkSmartPointer<vtkStructuredPoints>::New(); // the container which stores all fields
54  vtkSmartPointer<vtkMultiBlockDataSet> multiBDS =
55  vtkSmartPointer<vtkMultiBlockDataSet>::New ();
56  vtkSmartPointer<vtkXMLMultiBlockDataReader> reader =
57  vtkSmartPointer<vtkXMLMultiBlockDataReader>::New();
58  vtkSmartPointer<vtkXMLMultiBlockDataWriter> writer =
59  vtkSmartPointer<vtkXMLMultiBlockDataWriter>::New();
60 
61  densityfield->SetNumberOfComponents(1);
62  densityfield->SetName("density");
63  velocityfield->SetNumberOfComponents(3);
64  velocityfield->SetName("velocity");
65 
66  unsigned int xLatticeSize = (unsigned int)(this->physicalNodes_[0].size()/this->cdecoarseningfactor_);
67  unsigned int yLatticeSize = (unsigned int)(this->physicalNodes_.size()/this->cdecoarseningfactor_);
68 
69  uniformGrid->SetDimensions(xLatticeSize,yLatticeSize,1);
70  uniformGrid->SetOrigin(0,0,0);
71  uniformGrid->SetSpacing(this->cdecoarseningfactor_,
72  this->cdecoarseningfactor_,
73  this->cdecoarseningfactor_); // adjust spacing when coarsening
74 
75  for (auto ity=0; ity<this->physicalNodes_.size(); ity+=this->cdecoarseningfactor_) { // loop y
76  for (auto itx=0; itx<this->physicalNodes_[0].size(); itx+=this->cdecoarseningfactor_) { // loop x
77  densityfield->InsertNextValue(this->physicalNodes_[ity][itx]->getFluidSolver().getRho());
78  velocityfield->InsertNextTuple3(this->physicalNodes_[ity][itx]->getFluidSolver().getVelocity().x,
79  this->physicalNodes_[ity][itx]->getFluidSolver().getVelocity().y,
80  0.0);
81  }
82  }
83 
84  uniformGrid->GetPointData()->AddArray(densityfield);
85  uniformGrid->GetPointData()->AddArray(velocityfield);
86 
87  reader->SetFileName(filename.str().c_str());
88  struct stat buffer;
89  if (stat(filename.str().c_str(), &buffer) == 0) { // check if file already exists
90  reader->Update(); // if yes, load its content
91  multiBDS->ShallowCopy(reader->GetOutput());
92  }
93 
94  multiBDS->SetBlock(multiBDS->GetNumberOfBlocks(),uniformGrid);
95 
96  writer->SetFileName(filename.str().c_str());
97 #if VTK_MAJOR_VERSION <= 5
98  writer->SetInput(multiBDS);
99 #else
100  writer->SetInputData(multiBDS);
101 #endif
102  //writer->SetCompressorTypeToNone();
103  writer->SetDataModeToAscii();
104  //writer->SetDataModeToBinary();
105  writer->Write();
106 }
107 }
108 } // end namespace
const std::string filename_
filename_ Stores the filename of this functor
virtual void operator()(unsigned int time) const
operator() Writes the report