VertexMeshWriter.cpp

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 #include "VertexMeshWriter.hpp"
00030 #include "Version.hpp"
00031 
00035 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00036 struct MeshWriterIterators
00037 {
00039     typename AbstractMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator*  pNodeIter;
00041     typename VertexMesh<ELEMENT_DIM,SPACE_DIM>::VertexElementIterator* pElemIter;
00042 };
00043 
00045 // Implementation
00047 
00048 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00049 VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::VertexMeshWriter(const std::string& rDirectory,
00050                                                            const std::string& rBaseName,
00051                                                            const bool clearOutputDir)
00052     : AbstractMeshWriter<ELEMENT_DIM, SPACE_DIM>(rDirectory, rBaseName, clearOutputDir),
00053       mpMesh(NULL),
00054       mpIters(new MeshWriterIterators<ELEMENT_DIM,SPACE_DIM>),
00055       mpNodeMap(NULL),
00056       mNodeMapCurrentIndex(0)
00057 {
00058     mpIters->pNodeIter = NULL;
00059     mpIters->pElemIter = NULL;
00060 
00061 #ifdef CHASTE_VTK
00062      // Dubious, since we shouldn't yet know what any details of the mesh are.
00063      mpVtkUnstructedMesh = vtkUnstructuredGrid::New();
00064 #endif //CHASTE_VTK
00065 }
00066 
00067 
00068 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00069 VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::~VertexMeshWriter()
00070 {
00071     if (mpIters->pNodeIter)
00072     {
00073         delete mpIters->pNodeIter;
00074         delete mpIters->pElemIter;
00075     }
00076 
00077     delete mpIters;
00078 
00079     if (mpNodeMap)
00080     {
00081         delete mpNodeMap;
00082     }
00083 
00084 #ifdef CHASTE_VTK
00085      // Dubious, since we shouldn't yet know what any details of the mesh are.
00086      mpVtkUnstructedMesh->Delete(); // Reference counted
00087 #endif //CHASTE_VTK
00088 }
00089 
00090 
00091 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00092 std::vector<double> VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::GetNextNode()
00093 {
00094     if (mpMesh)
00095     {
00096         std::vector<double> coords(SPACE_DIM);
00097 
00098         assert(this->mNumNodes==mpMesh->GetNumNodes());
00099 
00100         // get the node coords using the node iterator (so to skip deleted nodes etc)
00101         for (unsigned j=0; j<SPACE_DIM; j++)
00102         {
00103             coords[j] = (*(mpIters->pNodeIter))->GetPoint()[j];
00104         }
00105 
00106         ++(*(mpIters->pNodeIter));
00107 
00108         return coords;
00109     }
00110     else
00111     {
00112         return AbstractMeshWriter<ELEMENT_DIM,SPACE_DIM>::GetNextNode();
00113     }
00114 }
00115 
00116 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00117 ElementData VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::GetNextElement()
00118 {
00119     if (mpMesh)
00120     {
00121         assert(this->mNumElements==mpMesh->GetNumElements());
00122 
00123         ElementData elem_data;
00124         elem_data.NodeIndices.resize((*(mpIters->pElemIter))->GetNumNodes());
00125         for (unsigned j=0; j<elem_data.NodeIndices.size(); j++)
00126         {
00127             unsigned old_index = (*(mpIters->pElemIter))->GetNodeGlobalIndex(j);
00128             elem_data.NodeIndices[j] = mpMesh->IsMeshChanging() ? mpNodeMap->GetNewIndex(old_index) : old_index;
00129         }
00130 // \todo: set attribute (#1076)
00131 
00132         ++(*(mpIters->pElemIter));
00133 
00134         return elem_data;
00135     }
00136     else
00137     {
00138         return AbstractMeshWriter<ELEMENT_DIM,SPACE_DIM>::GetNextElement();
00139     }
00140 }
00141 
00142 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00143 void VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::WriteVtkUsingMesh(VertexMesh<ELEMENT_DIM, SPACE_DIM>& rMesh, std::string stamp)
00144 {
00145 #ifdef CHASTE_VTK
00146     //Make the Vtk mesh
00147     assert(SPACE_DIM==3 || SPACE_DIM == 2);
00148     vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
00149     p_pts->GetData()->SetName("Vertex positions");
00150     for (unsigned node_num=0; node_num<rMesh.GetNumNodes(); node_num++)
00151     {
00152         c_vector<double, SPACE_DIM> position = rMesh.GetNode(node_num)->rGetLocation();
00153         if (SPACE_DIM==2)
00154         {
00155             p_pts->InsertPoint(node_num, position[0], position[1], 0.0);
00156         }
00157         else
00158         {
00159             p_pts->InsertPoint(node_num, position[0], position[1], position[2]);
00160         }
00161     }
00162 
00163     mpVtkUnstructedMesh->SetPoints(p_pts);
00164     p_pts->Delete(); //Reference counted
00165     for (typename VertexMesh<ELEMENT_DIM,SPACE_DIM>::VertexElementIterator iter = rMesh.GetElementIteratorBegin();
00166              iter != rMesh.GetElementIteratorEnd();
00167              ++iter)
00168     {
00169         vtkCell* p_cell;
00170         if (SPACE_DIM == 2)
00171         {
00172             p_cell = vtkPolygon::New();
00173         }
00174         else
00175         {
00176             p_cell = vtkConvexPointSet::New();
00177         }
00178         vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00179         p_cell_id_list->SetNumberOfIds(iter->GetNumNodes());
00180         for (unsigned j = 0; j < iter->GetNumNodes(); ++j)
00181         {
00182             p_cell_id_list->SetId(j, iter->GetNodeGlobalIndex(j));
00183         }
00184         mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00185         p_cell->Delete(); //Reference counted
00186     }
00187 
00188     //Vtk mesh is now made
00189     assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00190     vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
00191     p_writer->SetInput(mpVtkUnstructedMesh);
00192     //Uninitialised stuff arises (see #1079), but you can remove
00193     //valgrind problems by removing compression:
00194     // **** REMOVE WITH CAUTION *****
00195     p_writer->SetCompressor(NULL);
00196     // **** REMOVE WITH CAUTION *****
00197 
00198     std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName;
00199     if (stamp != "")
00200     {
00201         vtk_file_name += "_" + stamp;
00202     }
00203     vtk_file_name += ".vtu";
00204 
00205     p_writer->SetFileName(vtk_file_name.c_str());
00206     //p_writer->PrintSelf(std::cout, vtkIndent());
00207     p_writer->Write();
00208     p_writer->Delete(); //Reference counted
00209 #endif //CHASTE_VTK
00210 }
00211 
00212 
00213 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00214 void VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
00215 {
00216 #ifdef CHASTE_VTK
00217     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00218     p_scalars->SetName(dataName.c_str());
00219     for (unsigned i=0; i<dataPayload.size(); i++)
00220     {
00221         p_scalars->InsertNextValue(dataPayload[i]);
00222     }
00223 
00224     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00225     p_cell_data->AddArray(p_scalars);
00226     p_scalars->Delete(); //Reference counted
00227 #endif //CHASTE_VTK
00228 }
00229 
00230 
00231 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00232 void VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
00233 {
00234 #ifdef CHASTE_VTK
00235     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00236     p_scalars->SetName(dataName.c_str());
00237     for (unsigned i=0; i<dataPayload.size(); i++)
00238     {
00239         p_scalars->InsertNextValue(dataPayload[i]);
00240     }
00241 
00242     vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
00243     p_point_data->AddArray(p_scalars);
00244     p_scalars->Delete(); //Reference counted
00245 #endif //CHASTE_VTK
00246 }
00248 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00249 void VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::WriteFilesUsingMesh(VertexMesh<ELEMENT_DIM,SPACE_DIM>& rMesh)
00250 {
00251 
00252     this->mpMeshReader = NULL;
00253     mpMesh = &rMesh;
00254 
00255     this->mNumNodes = mpMesh->GetNumNodes();
00256     this->mNumElements = mpMesh->GetNumElements();
00257 
00258     typedef typename AbstractMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator NodeIterType;
00259     mpIters->pNodeIter = new NodeIterType(mpMesh->GetNodeIteratorBegin());
00260 
00261     typedef typename VertexMesh<ELEMENT_DIM,SPACE_DIM>::VertexElementIterator ElemIterType;
00262     mpIters->pElemIter = new ElemIterType(mpMesh->GetElementIteratorBegin());
00263 
00264     // Set up node map if we might have deleted nodes
00265     mNodeMapCurrentIndex = 0;
00266     if (mpMesh->IsMeshChanging())
00267     {
00268         mpNodeMap = new NodeMap(mpMesh->GetNumAllNodes());
00269         for (NodeIterType it = mpMesh->GetNodeIteratorBegin(); it != mpMesh->GetNodeIteratorEnd(); ++it)
00270         {
00271             mpNodeMap->SetNewIndex(it->GetIndex(), mNodeMapCurrentIndex++);
00272         }
00273     }
00274     WriteFiles();
00275 }
00276 
00277 
00278 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00279 void VertexMeshWriter<ELEMENT_DIM, SPACE_DIM>::WriteFiles()
00280 {
00281     std::string comment = "# " + ChasteBuildInfo::GetProvenanceString();
00282 
00283     // Write node file
00284     std::string node_file_name = this->mBaseName + ".node";
00285     out_stream p_node_file = this->mpOutputFileHandler->OpenOutputFile(node_file_name);
00286 
00287     // Write the node header
00288     unsigned num_attr = 0;
00289     unsigned max_bdy_marker = 0;
00290     unsigned num_nodes = this->GetNumNodes();
00291 
00292     *p_node_file << num_nodes << "\t";
00293     *p_node_file << SPACE_DIM << "\t";
00294     *p_node_file << num_attr << "\t";
00295     *p_node_file << max_bdy_marker << "\n";
00296     *p_node_file << std::setprecision(6);
00297 
00298     // Write each node's data
00299     unsigned default_marker = 0;
00300     for (unsigned item_num=0; item_num<num_nodes; item_num++)
00301     {
00302         std::vector<double> current_item = this->GetNextNode();
00303         *p_node_file << item_num;
00304         for (unsigned i=0; i<SPACE_DIM; i++)
00305         {
00306             *p_node_file << "\t" << current_item[i];
00307         }
00308         *p_node_file << "\t" << default_marker << "\n";
00309 
00310     }
00311     *p_node_file << comment << "\n";
00312     p_node_file->close();
00313 
00314     // Write element file
00315     std::string element_file_name = this->mBaseName + ".cell";
00316     out_stream p_element_file = this->mpOutputFileHandler->OpenOutputFile(element_file_name);
00317 
00318     // Write the element header
00319     unsigned num_elements = this->GetNumElements();
00320 
00321     *p_element_file << num_elements << "\t";
00322     *p_element_file << num_attr << "\n";
00323 
00324     // Write each element's data
00326     for (unsigned item_num=0; item_num<num_elements; item_num++)
00327     {
00328         std::vector<unsigned> current_item = this->GetNextElement().NodeIndices;
00329 
00330         *p_element_file << item_num <<  "\t" << current_item.size();
00331         for (unsigned i=0; i<current_item.size(); i++)
00332         {
00333             *p_element_file << "\t" << current_item[i];
00334         }
00335         *p_element_file << "\n";
00336     }
00337     *p_element_file << comment << "\n";
00338     p_element_file->close();
00339 }
00340 
00342 // Explicit instantiation
00344 
00345 template class VertexMeshWriter<1,1>;
00346 template class VertexMeshWriter<1,2>;
00347 template class VertexMeshWriter<1,3>;
00348 template class VertexMeshWriter<2,2>;
00349 template class VertexMeshWriter<2,3>;
00350 template class VertexMeshWriter<3,3>;

Generated by  doxygen 1.6.2