VtkWriter.hpp
00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef VTKWRITER_HPP_
00030 #define VTKWRITER_HPP_
00031
00032 #ifdef CHASTE_VTK
00033
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
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
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
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();
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();
00121 if (SPACE_DIM==2)
00122 {
00123 current_item.push_back(0.0);
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();
00131 for (unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
00132 {
00133 std::vector<unsigned> current_element = this->GetNextElement().NodeIndices;
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();
00151 }
00152 }
00153
00154 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00155 void VtkWriter<ELEMENT_DIM,SPACE_DIM>::WriteFiles()
00156 {
00157
00158 {
00159 MakeVtkMesh();
00160 assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00161 vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
00162 p_writer->SetInput(mpVtkUnstructedMesh);
00163
00164
00165
00166 p_writer->SetCompressor(NULL);
00167
00168 std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
00169 p_writer->SetFileName(vtk_file_name.c_str());
00170
00171 p_writer->Write();
00172 p_writer->Delete();
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();
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();
00213
00214 }
00215 #endif //CHASTE_VTK
00216
00217 #endif