Hdf5ToVtkConverter.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
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 "UblasCustomFunctions.hpp"
00030 #include "Hdf5ToVtkConverter.hpp"
00031 #include "PetscTools.hpp"
00032 #include "Exception.hpp"
00033 #include "ReplicatableVector.hpp"
00034 #include "DistributedVector.hpp"
00035 #include "DistributedVectorFactory.hpp"
00036 #include "VtkMeshWriter.hpp"
00037 #include "GenericMeshReader.hpp"
00038 #include "DistributedTetrahedralMesh.hpp"
00039 #include "Warnings.hpp"
00040 
00041 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00042 Hdf5ToVtkConverter<ELEMENT_DIM, SPACE_DIM>::Hdf5ToVtkConverter(std::string inputDirectory,
00043                                                                std::string fileBaseName,
00044                                                                AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00045                                                                bool parallelVtk,
00046                                                                bool usingOriginalNodeOrdering)
00047     : AbstractHdf5Converter<ELEMENT_DIM,SPACE_DIM>(inputDirectory, fileBaseName, pMesh, "vtk_output")
00048 {
00049 #ifdef CHASTE_VTK // Requires "sudo aptitude install libvtk5-dev" or similar
00050 
00051     // Write mesh in a suitable form for VTK
00052     std::string output_directory = inputDirectory + "/" + this->mRelativeSubdirectory;
00053     VtkMeshWriter<ELEMENT_DIM,SPACE_DIM> vtk_writer(output_directory, fileBaseName, false);
00054 
00055     DistributedVectorFactory* p_factory = pMesh->GetDistributedVectorFactory();
00056 
00057     // Make sure that we are never trying to write from an incomplete data HDF5 file
00058     assert(this->mpReader->GetNumberOfRows() == pMesh->GetNumNodes());
00059 
00060     DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* p_distributed_mesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>*>(pMesh);
00061 
00062     unsigned num_nodes = pMesh->GetNumNodes();
00063     if (parallelVtk)
00064     {
00065         // If it's not a distributed mesh, then we might want to give a warning and back-off
00066         if (p_distributed_mesh == NULL)
00067         {
00068             WARNING("Can only write parallel VTK from a DistributedTetrahedralMesh - writing sequential VTK instead");
00069             parallelVtk = false;
00070         }
00071 
00072         // If the node ordering flag is set, then we can't do this
00073         if (usingOriginalNodeOrdering)
00074         {
00075             WARNING("Can't write parallel VTK (pvtu) files with original ordering - writing sequential VTK instead");
00076             parallelVtk = false;
00077         }
00078 
00079         // Are we now committed to writing .pvtu?
00080         if (parallelVtk)
00081         {
00082            vtk_writer.SetParallelFiles(*pMesh);
00083            num_nodes = p_distributed_mesh->GetNumLocalNodes();
00084         }
00085     }
00086 
00087     Vec data = p_factory->CreateVec();
00088 
00089     unsigned num_timesteps = this->mpReader->GetUnlimitedDimensionValues().size();
00090 
00091     // Loop over time steps
00092     for (unsigned time_step=0; time_step<num_timesteps; time_step++)
00093     {
00094         // Loop over variables
00095         for (unsigned variable=0; variable<this->mNumVariables; variable++)
00096         {
00097             std::string variable_name = this->mpReader->GetVariableNames()[variable];
00098 
00099             // Gets variable at this time step from HDF5 archive
00100             this->mpReader->GetVariableOverNodes(data, variable_name, time_step);
00101 
00102             std::vector<double> data_for_vtk;
00103             data_for_vtk.resize(num_nodes);
00104             std::ostringstream variable_point_data_name;
00105             variable_point_data_name << variable_name << "_" << std::setw(6) << std::setfill('0') << time_step;
00106 
00107             if (parallelVtk)
00108             {
00109                 // Parallel VTU files
00110                 double *p_data;
00111                 VecGetArray(data, &p_data);
00112                 for (unsigned index=0; index<num_nodes; index++)
00113                 {
00114                     data_for_vtk[index]  = p_data[index];
00115                 }
00116                 VecRestoreArray(data, &p_data);
00117             }
00118             else
00119             {
00120                 // One VTU file
00121                 ReplicatableVector repl_data(data);
00122                 for (unsigned index=0; index<num_nodes; index++)
00123                 {
00124                     data_for_vtk[index] = repl_data[index];
00125                 }
00126             }
00127             // Add this variable into the node "point" data
00128             vtk_writer.AddPointData(variable_point_data_name.str(), data_for_vtk);
00129         }
00130     }
00131 
00132     // Tidy up
00133     VecDestroy(data);
00134 
00135     // Normally the in-memory mesh is converted
00136     if (!usingOriginalNodeOrdering)
00137     {
00138         vtk_writer.WriteFilesUsingMesh(*(this->mpMesh));
00139     }
00140     else
00141     {
00142         // In this case we expect the mesh to have been read in from file
00144         // Note that the next line will throw if the mesh has not been read from file
00145         std::string original_file = this->mpMesh->GetMeshFileBaseName();
00146         std::auto_ptr<AbstractMeshReader<ELEMENT_DIM, SPACE_DIM> > p_original_mesh_reader
00147             = GenericMeshReader<ELEMENT_DIM, SPACE_DIM>(original_file);
00148         vtk_writer.WriteFilesUsingMeshReader(*p_original_mesh_reader);
00149     }
00150 #endif //CHASTE_VTK
00151 }
00152 
00154 // Explicit instantiation
00156 
00157 template class Hdf5ToVtkConverter<1,1>;
00158 template class Hdf5ToVtkConverter<1,2>;
00159 template class Hdf5ToVtkConverter<2,2>;
00160 template class Hdf5ToVtkConverter<1,3>;
00161 template class Hdf5ToVtkConverter<2,3>;
00162 template class Hdf5ToVtkConverter<3,3>;
Generated on Thu Dec 22 13:00:17 2011 for Chaste by  doxygen 1.6.3