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

Generated on Mon Apr 18 11:35:28 2011 for Chaste by  doxygen 1.5.5