LBIBCell
 All Classes Functions Variables Friends Pages
vtkCDEReporter.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/vtkCDEReporter.hpp>
23 #include <UtilLib/include/Exception.hpp>
24 #include <LbmLib/include/nodes/PhysicalNode.hpp>
25 #include <LbmLib/include/solver/CDESolver/CDEAbstractSolver.hpp>
26 
27 #include <vtkSmartPointer.h>
28 #include <vtkDoubleArray.h>
29 #include <vtkPointData.h>
30 
31 #include <vtkMultiBlockDataSet.h>
32 #include <vtkXMLMultiBlockDataWriter.h>
33 #include <vtkXMLMultiBlockDataReader.h>
34 #include <vtkStructuredPoints.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 vtkCDEReporter::operator()(unsigned int time) const {
44  std::stringstream filename;
45  filename << filename_ << "_" << time << ".vtm";
46 
47  std::vector<solver::CDEAbstractSolver*> cdesolvers;
48  std::map<std::string,vtkSmartPointer<vtkDoubleArray> > cdesolvermap;
49  vtkSmartPointer<vtkMultiBlockDataSet> multiBDS =
50  vtkSmartPointer<vtkMultiBlockDataSet>::New ();
51  vtkSmartPointer<vtkXMLMultiBlockDataReader> reader =
52  vtkSmartPointer<vtkXMLMultiBlockDataReader>::New();
53  vtkSmartPointer<vtkXMLMultiBlockDataWriter> writer =
54  vtkSmartPointer<vtkXMLMultiBlockDataWriter>::New();
55  vtkSmartPointer<vtkStructuredPoints> uniformGrid =
56  vtkSmartPointer<vtkStructuredPoints>::New(); // the container which stores all fields
57 
58  cdesolvers = this->physicalNodes_[0][0]->getCDESolvers();
59  unsigned int numcdesolvers = cdesolvers.size(); // number of cde solvers. consistency check later.
60  for (auto k=0; k<numcdesolvers; ++k) {
61  cdesolvermap[cdesolvers[k]->getName()] = vtkSmartPointer<vtkDoubleArray>::New();
62  cdesolvermap[cdesolvers[k]->getName()]->SetNumberOfComponents(1);
63  cdesolvermap[cdesolvers[k]->getName()]->SetName(cdesolvers[k]->getName().c_str());
64  }
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 
70  uniformGrid->SetDimensions(xLatticeSize,yLatticeSize,1);
71  uniformGrid->SetOrigin(0,0,0);
72  uniformGrid->SetSpacing(this->cdecoarseningfactor_,
73  this->cdecoarseningfactor_,
74  this->cdecoarseningfactor_); // adjust spacing when coarsening
75 
76  for (auto ity=0; ity<this->physicalNodes_.size(); ity+=this->cdecoarseningfactor_) { // loop y
77  for (auto itx=0; itx<this->physicalNodes_[0].size(); itx+=this->cdecoarseningfactor_) { // loop x
78  cdesolvers =this->physicalNodes_[ity][itx]->getCDESolvers(); // get all the solvers from node [y][x]
79  for (auto k=0; k<numcdesolvers; ++k) { // loop all cdesolvers
80  cdesolvermap[cdesolvers[k]->getName()]->InsertNextValue(cdesolvers[k]->getC()); // add the concentration
81  }
82  }
83  }
84 
85  for (auto it : cdesolvermap) { // add all the fields to the uniform grid:
86  uniformGrid->GetPointData()->AddArray(it.second);
87  }
88 
89  reader->SetFileName(filename.str().c_str());
90  struct stat buffer;
91  if (stat(filename.str().c_str(), &buffer) == 0) { // check if already existing
92  reader->Update(); // if yes, load its content
93  multiBDS->ShallowCopy(reader->GetOutput());
94  }
95 
96  multiBDS->SetBlock(multiBDS->GetNumberOfBlocks(),uniformGrid);
97 
98  writer->SetFileName(filename.str().c_str());
99 #if VTK_MAJOR_VERSION <= 5
100  writer->SetInput(multiBDS);
101 #else
102  writer->SetInputData(multiBDS);
103 #endif
104  //writer->SetCompressorTypeToNone();
105  writer->SetDataModeToAscii();
106  //writer->SetDataModeToBinary();
107  writer->Write();
108 }
109 }
110 } // end namespace
virtual void operator()(unsigned int time) const
operator() Writes the report
const std::string filename_
filename_ Stores the filename of this functor