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 VTKMESHWRITER_HPP_
00030 #define VTKMESHWRITER_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 VtkMeshWriter : public AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>
00056 {
00057
00058
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
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
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();
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();
00145 if (SPACE_DIM==2)
00146 {
00147 current_item.push_back(0.0);
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();
00155 for (unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
00156 {
00157 std::vector<unsigned> current_element = this->GetNextElement().NodeIndices;
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();
00175 }
00176 }
00177
00178 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00179 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::WriteFiles()
00180 {
00181
00182 {
00183 MakeVtkMesh();
00184 assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00185 vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
00186 p_writer->SetInput(mpVtkUnstructedMesh);
00187
00188
00189
00190 p_writer->SetCompressor(NULL);
00191
00192 std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
00193 p_writer->SetFileName(vtk_file_name.c_str());
00194
00195 p_writer->Write();
00196 p_writer->Delete();
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();
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();
00237
00238 }
00239 #endif //CHASTE_VTK
00240
00241 #endif