VtkMeshWriter.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 VTKMESHWRITER_HPP_
00030 #define VTKMESHWRITER_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 VtkMeshWriter : public AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>
00056 {
00057 
00058 //Requires  "sudo aptitude install libvtk5-dev" or similar
00059 
00060 private:
00067     vtkUnstructuredGrid* mpVtkUnstructedMesh;
00068 
00073     void MakeVtkMesh();
00074 public:
00075 
00083     VtkMeshWriter(const std::string& rDirectory, const std::string& rBaseName, const bool& rCleanDirectory=true);
00084 
00088     void WriteFiles();
00089 
00097     void AddCellData(std::string name, std::vector<double> data);
00105     void AddPointData(std::string name, std::vector<double> data);
00106 
00110     virtual ~VtkMeshWriter();
00111 };
00112 
00113 
00115 // Implementation
00117 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00118 VtkMeshWriter<ELEMENT_DIM, SPACE_DIM>::VtkMeshWriter(const std::string& rDirectory,
00119                      const std::string& rBaseName,
00120                      const bool& rCleanDirectory)
00121     : AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>(rDirectory, rBaseName, rCleanDirectory)
00122 {
00123     this->mIndexFromZero = true;
00124 
00125     // Dubious, since we shouldn't yet know what any details of the mesh are.
00126     mpVtkUnstructedMesh = vtkUnstructuredGrid::New();
00127 }
00128 
00129 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00130 VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::~VtkMeshWriter()
00131 {
00132     mpVtkUnstructedMesh->Delete(); // Reference counted
00133 }
00134 
00135 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00136 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::MakeVtkMesh()
00137 {
00138     assert(SPACE_DIM==3 || SPACE_DIM == 2);
00139     assert(SPACE_DIM==ELEMENT_DIM);
00140     vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
00141     p_pts->GetData()->SetName("Vertex positions");
00142     for (unsigned item_num=0; item_num<this->GetNumNodes(); item_num++)
00143     {
00144         std::vector<double> current_item = this->GetNextNode(); //this->mNodeData[item_num];
00145         if (SPACE_DIM==2)
00146         {
00147             current_item.push_back(0.0);//For z-coordinate
00148         }
00149         assert(current_item.size() == 3);
00150         p_pts->InsertPoint(item_num, current_item[0], current_item[1], current_item[2]);
00151     }
00152 
00153     mpVtkUnstructedMesh->SetPoints(p_pts);
00154     p_pts->Delete(); //Reference counted
00155     for (unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
00156     {
00157         std::vector<unsigned> current_element = this->GetNextElement().NodeIndices; // this->mElementData[item_num];
00158         assert(current_element.size() == ELEMENT_DIM + 1);
00159         vtkCell* p_cell=NULL;
00160         if (SPACE_DIM == 3)
00161         {
00162             p_cell = vtkTetra::New();
00163         }
00164         if (SPACE_DIM == 2)
00165         {
00166             p_cell = vtkTriangle::New();
00167         }
00168         vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00169         for (unsigned j = 0; j < ELEMENT_DIM+1; ++j)
00170         {
00171             p_cell_id_list->SetId(j, current_element[j]);
00172         }
00173         mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00174         p_cell->Delete(); //Reference counted
00175     }
00176 }
00177 
00178 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00179 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::WriteFiles()
00180 {
00181     // Using separate scope here to make sure file is properly closed before re-opening it to add provenance info.
00182     {
00183         MakeVtkMesh();
00184         assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00185         vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
00186         p_writer->SetInput(mpVtkUnstructedMesh);
00187         //Uninitialised stuff arises (see #1079), but you can remove
00188         //valgrind problems by removing compression:
00189         // **** REMOVE WITH CAUTION *****
00190         p_writer->SetCompressor(NULL);
00191         // **** REMOVE WITH CAUTION *****
00192         std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
00193         p_writer->SetFileName(vtk_file_name.c_str());
00194         //p_writer->PrintSelf(std::cout, vtkIndent());
00195         p_writer->Write();
00196         p_writer->Delete(); //Reference counted
00197     }    
00198 
00199     std::string comment = "<!-- " + ChasteBuildInfo::GetProvenanceString() + "-->";
00200     
00201     std::string node_file_name = this->mBaseName + ".vtu";
00202     out_stream p_node_file = this->mpOutputFileHandler->OpenOutputFile(node_file_name, std::ios::out | std::ios::app);
00203 
00204     *p_node_file << "\n" << comment << "\n";    
00205     p_node_file->close();
00206     
00207 }
00208 
00209 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00210 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
00211 {
00212     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00213     p_scalars->SetName(dataName.c_str());
00214     for (unsigned i=0; i<dataPayload.size(); i++)
00215     {
00216         p_scalars->InsertNextValue(dataPayload[i]);
00217     }
00218 
00219     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00220     p_cell_data->AddArray(p_scalars);
00221     p_scalars->Delete(); //Reference counted
00222 }
00223 
00224 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00225 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
00226 {
00227     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00228     p_scalars->SetName(dataName.c_str());
00229     for (unsigned i=0; i<dataPayload.size(); i++)
00230     {
00231         p_scalars->InsertNextValue(dataPayload[i]);
00232     }
00233 
00234     vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
00235     p_point_data->AddArray(p_scalars);
00236     p_scalars->Delete(); //Reference counted
00237 
00238 }
00239 #endif //CHASTE_VTK
00240 
00241 #endif /*VTKMESHWRITER_HPP_*/

Generated on Mon Nov 1 12:35:24 2010 for Chaste by  doxygen 1.5.5