VtkWriter.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 
00053 template <unsigned DIM>
00054 class VtkWriter : public AbstractTetrahedralMeshWriter<DIM, DIM>
00055 {
00056 
00057 //Requires  "sudo aptitude install libvtk5-dev" or similar
00058 
00059 private:
00060     vtkUnstructuredGrid *mpVtkUnstructedMesh;
00061 
00062     void MakeVtkMesh();
00063 public:
00064 
00072     VtkWriter(const std::string& rDirectory, const std::string& rBaseName, const bool& rCleanDirectory=true);
00073 
00077     void WriteFiles();
00078 
00079     void AddCellData(std::string name, std::vector<double> data);
00080     void AddPointData(std::string name, std::vector<double> data);
00081 
00085     virtual ~VtkWriter();
00086 };
00087 
00088 
00090 // Implementation
00092 template <unsigned DIM>
00093 VtkWriter<DIM>::VtkWriter(const std::string& rDirectory,
00094                      const std::string& rBaseName,
00095                      const bool& rCleanDirectory)
00096     : AbstractTetrahedralMeshWriter<DIM, DIM>(rDirectory, rBaseName, rCleanDirectory)
00097 {
00098     this->mIndexFromZero = true;
00099 
00100     // Dubious, since we shouldn't yet know what any details of the mesh are.
00101     mpVtkUnstructedMesh = vtkUnstructuredGrid::New();
00102 }
00103 
00104 template <unsigned DIM>
00105 VtkWriter<DIM>::~VtkWriter()
00106 {
00107     mpVtkUnstructedMesh->Delete(); // Reference counted
00108 }
00109 
00110 template <unsigned DIM>
00111 void VtkWriter<DIM>::MakeVtkMesh()
00112 {
00113     assert(DIM==3 || DIM == 2);
00114     vtkPoints *p_pts = vtkPoints::New(VTK_DOUBLE);
00115     //p_pts->SetDataTypeToDouble();
00116     p_pts->GetData()->SetName("Vertex positions");
00117     for (unsigned item_num=0; item_num<this->GetNumNodes(); item_num++)
00118     {
00119         std::vector<double> current_item = this->mNodeData[item_num];
00120         if (DIM==2)
00121         {
00122             current_item.push_back(0.0);//For z-coordinate
00123         }
00124         assert(current_item.size() == 3);
00125         p_pts->InsertPoint(item_num, current_item[0], current_item[1], current_item[2]);
00126     }
00127 
00128     //mpVtkUnstructedMesh->Allocate(rMesh.GetNumNodes(), rMesh.GetNumNodes());
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->mElementData[item_num];
00134         assert(current_element.size() == DIM + 1);
00135         vtkCell *p_cell;
00136         if (DIM == 3)
00137         {
00138             p_cell = vtkTetra::New();
00139         }
00140         if (DIM == 2)
00141         {
00142             p_cell = vtkTriangle::New();
00143         }
00144         vtkIdList *p_cell_id_list = p_cell->GetPointIds();
00145         for (unsigned j = 0; j < 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 DIM>
00155 void VtkWriter<DIM>::WriteFiles()
00156 {
00157     MakeVtkMesh();
00158     assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00159     vtkXMLUnstructuredGridWriter *p_writer = vtkXMLUnstructuredGridWriter::New();
00160     p_writer->SetInput(mpVtkUnstructedMesh);
00161     //Uninitialised stuff arises (see #1079), but you can remove
00162     //valgrind problems by removing compression:
00163     // **** REMOVE WITH CAUTION *****
00164     p_writer->SetCompressor(NULL);
00165     // **** REMOVE WITH CAUTION *****
00166     std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
00167     p_writer->SetFileName(vtk_file_name.c_str());
00168     //p_writer->PrintSelf(std::cout, vtkIndent());
00169     p_writer->Write();
00170     p_writer->Delete(); //Reference counted
00171 }
00172 
00173 template <unsigned DIM>
00174 void VtkWriter<DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
00175 {
00176     vtkDoubleArray *p_scalars = vtkDoubleArray::New();
00177     p_scalars->SetName(dataName.c_str());
00178     for (unsigned i=0; i<dataPayload.size(); i++)
00179     {
00180         p_scalars->InsertNextValue(dataPayload[i]);
00181     }
00182 
00183     vtkCellData *p_cell_data = mpVtkUnstructedMesh->GetCellData();
00184     p_cell_data->AddArray(p_scalars);
00185     p_scalars->Delete(); //Reference counted
00186 }
00187 
00188 template <unsigned DIM>
00189 void VtkWriter<DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
00190 {
00191     vtkDoubleArray *p_scalars = vtkDoubleArray::New();
00192     p_scalars->SetName(dataName.c_str());
00193     for (unsigned i=0; i<dataPayload.size(); i++)
00194     {
00195         p_scalars->InsertNextValue(dataPayload[i]);
00196     }
00197 
00198     vtkPointData *p_point_data = mpVtkUnstructedMesh->GetPointData();
00199     p_point_data->AddArray(p_scalars);
00200     p_scalars->Delete(); //Reference counted
00201 
00202 }
00203 #endif //CHASTE_VTK
00204 
00205 #endif /*VTKWRITER_HPP_*/

Generated on Tue Aug 4 16:10:24 2009 for Chaste by  doxygen 1.5.5