VtkWriter.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #ifndef VTKWRITER_HPP_
00030 #define VTKWRITER_HPP_
00031 
00032 #ifdef CHASTE_VTK
00033 //Requires  "sudo aptitude install libvtk5-dev" or similar
00034 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the strstream deprecated warning for now (gcc4.3)
00035 #include <vtkDoubleArray.h>
00036 #include <vtkCellData.h>
00037 #include <vtkPointData.h>
00038 #include <vtkTetra.h>
00039 #include <vtkTriangle.h>
00040 #include <vtkUnstructuredGrid.h>
00041 #include <vtkUnstructuredGridWriter.h>
00042 #include <vtkXMLUnstructuredGridWriter.h>
00043 
00044 #include <vtkDataCompressor.h>
00045 #include "AbstractTetrahedralMeshWriter.hpp"
00046 #include "Version.hpp"
00047 
00054 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00055 class VtkWriter : public AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>
00056 {
00057 
00058 //Requires  "sudo aptitude install libvtk5-dev" or similar
00059 
00060 private:
00061     vtkUnstructuredGrid* mpVtkUnstructedMesh;
00062 
00063     void MakeVtkMesh();
00064 public:
00065 
00073     VtkWriter(const std::string& rDirectory, const std::string& rBaseName, const bool& rCleanDirectory=true);
00074 
00078     void WriteFiles();
00079 
00080     void AddCellData(std::string name, std::vector<double> data);
00081     void AddPointData(std::string name, std::vector<double> data);
00082 
00086     virtual ~VtkWriter();
00087 };
00088 
00089 
00091 // Implementation
00093 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00094 VtkWriter<ELEMENT_DIM, SPACE_DIM>::VtkWriter(const std::string& rDirectory,
00095                      const std::string& rBaseName,
00096                      const bool& rCleanDirectory)
00097     : AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>(rDirectory, rBaseName, rCleanDirectory)
00098 {
00099     this->mIndexFromZero = true;
00100 
00101     // Dubious, since we shouldn't yet know what any details of the mesh are.
00102     mpVtkUnstructedMesh = vtkUnstructuredGrid::New();
00103 }
00104 
00105 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00106 VtkWriter<ELEMENT_DIM,SPACE_DIM>::~VtkWriter()
00107 {
00108     mpVtkUnstructedMesh->Delete(); // Reference counted
00109 }
00110 
00111 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00112 void VtkWriter<ELEMENT_DIM,SPACE_DIM>::MakeVtkMesh()
00113 {
00114     assert(SPACE_DIM==3 || SPACE_DIM == 2);
00115     assert(SPACE_DIM==ELEMENT_DIM);
00116     vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
00117     p_pts->GetData()->SetName("Vertex positions");
00118     for (unsigned item_num=0; item_num<this->GetNumNodes(); item_num++)
00119     {
00120         std::vector<double> current_item = this->GetNextNode(); //this->mNodeData[item_num];
00121         if (SPACE_DIM==2)
00122         {
00123             current_item.push_back(0.0);//For z-coordinate
00124         }
00125         assert(current_item.size() == 3);
00126         p_pts->InsertPoint(item_num, current_item[0], current_item[1], current_item[2]);
00127     }
00128 
00129     mpVtkUnstructedMesh->SetPoints(p_pts);
00130     p_pts->Delete(); //Reference counted
00131     for (unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
00132     {
00133         std::vector<unsigned> current_element = this->GetNextElement().NodeIndices; // this->mElementData[item_num];
00134         assert(current_element.size() == ELEMENT_DIM + 1);
00135         vtkCell* p_cell=NULL;
00136         if (SPACE_DIM == 3)
00137         {
00138             p_cell = vtkTetra::New();
00139         }
00140         if (SPACE_DIM == 2)
00141         {
00142             p_cell = vtkTriangle::New();
00143         }
00144         vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00145         for (unsigned j = 0; j < ELEMENT_DIM+1; ++j)
00146         {
00147             p_cell_id_list->SetId(j, current_element[j]);
00148         }
00149         mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00150         p_cell->Delete(); //Reference counted
00151     }
00152 }
00153 
00154 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00155 void VtkWriter<ELEMENT_DIM,SPACE_DIM>::WriteFiles()
00156 {
00157     // Using separate scope here to make sure file is properly closed before re-opening it to add provenance info.
00158     {
00159         MakeVtkMesh();
00160         assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00161         vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
00162         p_writer->SetInput(mpVtkUnstructedMesh);
00163         //Uninitialised stuff arises (see #1079), but you can remove
00164         //valgrind problems by removing compression:
00165         // **** REMOVE WITH CAUTION *****
00166         p_writer->SetCompressor(NULL);
00167         // **** REMOVE WITH CAUTION *****
00168         std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
00169         p_writer->SetFileName(vtk_file_name.c_str());
00170         //p_writer->PrintSelf(std::cout, vtkIndent());
00171         p_writer->Write();
00172         p_writer->Delete(); //Reference counted
00173     }    
00174 
00175     std::string comment = "<!-- " + ChasteBuildInfo::GetProvenanceString() + "-->";
00176     
00177     std::string node_file_name = this->mBaseName + ".vtu";
00178     out_stream p_node_file = this->mpOutputFileHandler->OpenOutputFile(node_file_name, std::ios::out | std::ios::app);
00179 
00180     *p_node_file << "\n" << comment << "\n";    
00181     p_node_file->close();
00182     
00183 }
00184 
00185 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00186 void VtkWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
00187 {
00188     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00189     p_scalars->SetName(dataName.c_str());
00190     for (unsigned i=0; i<dataPayload.size(); i++)
00191     {
00192         p_scalars->InsertNextValue(dataPayload[i]);
00193     }
00194 
00195     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00196     p_cell_data->AddArray(p_scalars);
00197     p_scalars->Delete(); //Reference counted
00198 }
00199 
00200 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00201 void VtkWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
00202 {
00203     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00204     p_scalars->SetName(dataName.c_str());
00205     for (unsigned i=0; i<dataPayload.size(); i++)
00206     {
00207         p_scalars->InsertNextValue(dataPayload[i]);
00208     }
00209 
00210     vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
00211     p_point_data->AddArray(p_scalars);
00212     p_scalars->Delete(); //Reference counted
00213 
00214 }
00215 #endif //CHASTE_VTK
00216 
00217 #endif /*VTKWRITER_HPP_*/

Generated by  doxygen 1.6.2