Chaste Release::3.1
Hdf5ToVtkConverter.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include "UblasCustomFunctions.hpp"
00037 #include "Hdf5ToVtkConverter.hpp"
00038 #include "PetscTools.hpp"
00039 #include "Exception.hpp"
00040 #include "ReplicatableVector.hpp"
00041 #include "DistributedVector.hpp"
00042 #include "DistributedVectorFactory.hpp"
00043 #include "VtkMeshWriter.hpp"
00044 #include "GenericMeshReader.hpp"
00045 #include "DistributedTetrahedralMesh.hpp"
00046 #include "Warnings.hpp"
00047 
00048 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00049 Hdf5ToVtkConverter<ELEMENT_DIM, SPACE_DIM>::Hdf5ToVtkConverter(std::string inputDirectory,
00050                                                                std::string fileBaseName,
00051                                                                AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00052                                                                bool parallelVtk,
00053                                                                bool usingOriginalNodeOrdering)
00054     : AbstractHdf5Converter<ELEMENT_DIM,SPACE_DIM>(inputDirectory, fileBaseName, pMesh, "vtk_output")
00055 {
00056 #ifdef CHASTE_VTK // Requires "sudo aptitude install libvtk5-dev" or similar
00057 
00058     // Write mesh in a suitable form for VTK
00059     std::string output_directory = inputDirectory + "/" + this->mRelativeSubdirectory;
00060     VtkMeshWriter<ELEMENT_DIM,SPACE_DIM> vtk_writer(output_directory, fileBaseName, false);
00061 
00062     DistributedVectorFactory* p_factory = pMesh->GetDistributedVectorFactory();
00063 
00064     // Make sure that we are never trying to write from an incomplete data HDF5 file
00065     assert(this->mpReader->GetNumberOfRows() == pMesh->GetNumNodes());
00066 
00067     DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* p_distributed_mesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>*>(pMesh);
00068 
00069     unsigned num_nodes = pMesh->GetNumNodes();
00070     if (parallelVtk)
00071     {
00072         // If it's not a distributed mesh, then we might want to give a warning and back-off
00073         if (p_distributed_mesh == NULL)
00074         {
00075             WARNING("Can only write parallel VTK from a DistributedTetrahedralMesh - writing sequential VTK instead");
00076             parallelVtk = false;
00077         }
00078 
00079         // If the node ordering flag is set, then we can't do this
00080         if (usingOriginalNodeOrdering)
00081         {
00082             WARNING("Can't write parallel VTK (pvtu) files with original ordering - writing sequential VTK instead");
00083             parallelVtk = false;
00084         }
00085 
00086         // Are we now committed to writing .pvtu?
00087         if (parallelVtk)
00088         {
00089            vtk_writer.SetParallelFiles(*pMesh);
00090            num_nodes = p_distributed_mesh->GetNumLocalNodes();
00091         }
00092     }
00093 
00094     Vec data = p_factory->CreateVec();
00095 
00096     unsigned num_timesteps = this->mpReader->GetUnlimitedDimensionValues().size();
00097 
00098     // Loop over time steps
00099     for (unsigned time_step=0; time_step<num_timesteps; time_step++)
00100     {
00101         // Loop over variables
00102         for (unsigned variable=0; variable<this->mNumVariables; variable++)
00103         {
00104             std::string variable_name = this->mpReader->GetVariableNames()[variable];
00105 
00106             // Gets variable at this time step from HDF5 archive
00107             this->mpReader->GetVariableOverNodes(data, variable_name, time_step);
00108 
00109             std::vector<double> data_for_vtk;
00110             data_for_vtk.resize(num_nodes);
00111             std::ostringstream variable_point_data_name;
00112             variable_point_data_name << variable_name << "_" << std::setw(6) << std::setfill('0') << time_step;
00113 
00114             if (parallelVtk)
00115             {
00116                 // Parallel VTU files
00117                 double *p_data;
00118                 VecGetArray(data, &p_data);
00119                 for (unsigned index=0; index<num_nodes; index++)
00120                 {
00121                     data_for_vtk[index]  = p_data[index];
00122                 }
00123                 VecRestoreArray(data, &p_data);
00124             }
00125             else
00126             {
00127                 // One VTU file
00128                 ReplicatableVector repl_data(data);
00129                 for (unsigned index=0; index<num_nodes; index++)
00130                 {
00131                     data_for_vtk[index] = repl_data[index];
00132                 }
00133             }
00134             // Add this variable into the node "point" data
00135             vtk_writer.AddPointData(variable_point_data_name.str(), data_for_vtk);
00136         }
00137     }
00138 
00139     // Tidy up
00140     PetscTools::Destroy(data);
00141 
00142     // Normally the in-memory mesh is converted
00143     if (!usingOriginalNodeOrdering)
00144     {
00145         vtk_writer.WriteFilesUsingMesh(*(this->mpMesh));
00146     }
00147     else
00148     {
00149         // In this case we expect the mesh to have been read in from file
00151         // Note that the next line will throw if the mesh has not been read from file
00152         std::string original_file = this->mpMesh->GetMeshFileBaseName();
00153         std::auto_ptr<AbstractMeshReader<ELEMENT_DIM, SPACE_DIM> > p_original_mesh_reader
00154             = GenericMeshReader<ELEMENT_DIM, SPACE_DIM>(original_file);
00155         vtk_writer.WriteFilesUsingMeshReader(*p_original_mesh_reader);
00156     }
00157 #endif //CHASTE_VTK
00158 }
00159 
00161 // Explicit instantiation
00163 
00164 template class Hdf5ToVtkConverter<1,1>;
00165 template class Hdf5ToVtkConverter<1,2>;
00166 template class Hdf5ToVtkConverter<2,2>;
00167 template class Hdf5ToVtkConverter<1,3>;
00168 template class Hdf5ToVtkConverter<2,3>;
00169 template class Hdf5ToVtkConverter<3,3>;