Hdf5ToVtkConverter.cpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, 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(const FileFinder& rInputDirectory,
00050                                                                const std::string& rFileBaseName,
00051                                                                AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00052                                                                bool parallelVtk,
00053                                                                bool usingOriginalNodeOrdering)
00054     : AbstractHdf5Converter<ELEMENT_DIM,SPACE_DIM>(rInputDirectory, rFileBaseName, pMesh, "vtk_output",0u)
00055 {
00056 #ifdef CHASTE_VTK // Requires "sudo aptitude install libvtk5-dev" or similar
00057 
00058     // Write mesh in a suitable form for VTK
00059     FileFinder test_output("", RelativeTo::ChasteTestOutput);
00060     std::string output_directory = rInputDirectory.GetRelativePath(test_output) + "/" + this->mRelativeSubdirectory;
00061 
00062     VtkMeshWriter<ELEMENT_DIM,SPACE_DIM> vtk_writer(output_directory, rFileBaseName, false);
00063 
00064     DistributedVectorFactory* p_factory = pMesh->GetDistributedVectorFactory();
00065 
00066     // Make sure that we are never trying to write from an incomplete data HDF5 file
00067     assert(this->mpReader->GetNumberOfRows() == pMesh->GetNumNodes());
00068 
00069     DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* p_distributed_mesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>*>(pMesh);
00070 
00071     unsigned num_nodes = pMesh->GetNumNodes();
00072     if (parallelVtk)
00073     {
00074         // If it's not a distributed mesh, then we might want to give a warning and back-off
00075         if (p_distributed_mesh == NULL)
00076         {
00077             WARNING("Can only write parallel VTK from a DistributedTetrahedralMesh - writing sequential VTK instead");
00078             parallelVtk = false;
00079         }
00080 
00081         // If the node ordering flag is set, then we can't do this
00082         if (usingOriginalNodeOrdering)
00083         {
00084             WARNING("Can't write parallel VTK (pvtu) files with original ordering - writing sequential VTK instead");
00085             parallelVtk = false;
00086         }
00087 
00088         // Are we now committed to writing .pvtu?
00089         if (parallelVtk)
00090         {
00091            vtk_writer.SetParallelFiles(*pMesh);
00092            num_nodes = p_distributed_mesh->GetNumLocalNodes();
00093         }
00094     }
00095 
00096     Vec data = p_factory->CreateVec();
00097 
00098     do // Loop over datasets via MoveOntoNextDataset method in the abstract class
00099     {
00100         // Make sure that we are never trying to write from an incomplete HDF5 dataset.
00101         assert(this->mpReader->GetNumberOfRows() == pMesh->GetNumNodes());
00102 
00103         unsigned num_timesteps = this->mpReader->GetUnlimitedDimensionValues().size();
00104 
00105         // Loop over time steps
00106         for (unsigned time_step=0; time_step<num_timesteps; time_step++)
00107         {
00108             // Loop over variables
00109             for (unsigned variable=0; variable<this->mNumVariables; variable++)
00110             {
00111                 std::string variable_name = this->mpReader->GetVariableNames()[variable];
00112 
00113                 // Gets variable at this time step from HDF5 archive
00114                 this->mpReader->GetVariableOverNodes(data, variable_name, time_step);
00115 
00116                 std::vector<double> data_for_vtk;
00117                 data_for_vtk.resize(num_nodes);
00118                 std::ostringstream variable_point_data_name;
00119                 variable_point_data_name << variable_name << "_" << std::setw(6) << std::setfill('0') << time_step;
00120 
00121                 if (parallelVtk)
00122                 {
00123                     // Parallel VTU files
00124                     double *p_data;
00125                     VecGetArray(data, &p_data);
00126                     for (unsigned index=0; index<num_nodes; index++)
00127                     {
00128                         data_for_vtk[index]  = p_data[index];
00129                     }
00130                     VecRestoreArray(data, &p_data);
00131                 }
00132                 else
00133                 {
00134                     // One VTU file
00135                     ReplicatableVector repl_data(data);
00136                     for (unsigned index=0; index<num_nodes; index++)
00137                     {
00138                         data_for_vtk[index] = repl_data[index];
00139                     }
00140                 }
00141                 // Add this variable into the node "point" data
00142                 vtk_writer.AddPointData(variable_point_data_name.str(), data_for_vtk);
00143             }
00144         }
00145     }
00146     while ( this->MoveOntoNextDataset() );
00147 
00148     // Tidy up
00149     PetscTools::Destroy(data);
00150 
00151     // Normally the in-memory mesh is converted
00152     if (!usingOriginalNodeOrdering)
00153     {
00154         vtk_writer.WriteFilesUsingMesh(*(this->mpMesh));
00155     }
00156     else
00157     {
00158         // In this case we expect the mesh to have been read in from file
00160         // Note that the next line will throw if the mesh has not been read from file
00161         std::string original_file = this->mpMesh->GetMeshFileBaseName();
00162         std::auto_ptr<AbstractMeshReader<ELEMENT_DIM, SPACE_DIM> > p_original_mesh_reader
00163             = GenericMeshReader<ELEMENT_DIM, SPACE_DIM>(original_file);
00164         vtk_writer.WriteFilesUsingMeshReader(*p_original_mesh_reader);
00165     }
00166 #endif //CHASTE_VTK
00167 }
00168 
00170 // Explicit instantiation
00172 
00173 template class Hdf5ToVtkConverter<1,1>;
00174 template class Hdf5ToVtkConverter<1,2>;
00175 template class Hdf5ToVtkConverter<2,2>;
00176 template class Hdf5ToVtkConverter<1,3>;
00177 template class Hdf5ToVtkConverter<2,3>;
00178 template class Hdf5ToVtkConverter<3,3>;

Generated by  doxygen 1.6.2