Chaste  Release::2018.1
Hdf5ToVtkConverter.cpp
1 /*
2 
3 Copyright (c) 2005-2018, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "UblasCustomFunctions.hpp"
37 #include "Hdf5ToVtkConverter.hpp"
38 #include "PetscTools.hpp"
39 #include "Exception.hpp"
40 #include "ReplicatableVector.hpp"
41 #include "DistributedVector.hpp"
42 #include "DistributedVectorFactory.hpp"
43 #include "VtkMeshWriter.hpp"
44 #include "GenericMeshReader.hpp"
45 #include "DistributedTetrahedralMesh.hpp"
46 #include "Warnings.hpp"
47 
48 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
50  const std::string& rFileBaseName,
52  bool parallelVtk,
53  bool usingOriginalNodeOrdering)
54  : AbstractHdf5Converter<ELEMENT_DIM,SPACE_DIM>(rInputDirectory, rFileBaseName, pMesh, "vtk_output",0u)
55 {
56 #ifdef CHASTE_VTK // Requires "sudo aptitude install libvtk5-dev" or similar
57 
58  // Write mesh in a suitable form for VTK
60  std::string output_directory = rInputDirectory.GetRelativePath(test_output) + "/" + this->mRelativeSubdirectory;
61 
62  VtkMeshWriter<ELEMENT_DIM,SPACE_DIM> vtk_writer(output_directory, rFileBaseName, false);
63 
65 
66  // Make sure that we are never trying to write from an incomplete data HDF5 file
67  assert(this->mpReader->GetNumberOfRows() == pMesh->GetNumNodes());
68 
70 
71  unsigned num_nodes = pMesh->GetNumNodes();
72  if (parallelVtk)
73  {
74  // If it's not a distributed mesh, then we might want to give a warning and back-off
75  if (p_distributed_mesh == nullptr)
76  {
77  WARNING("Can only write parallel VTK from a DistributedTetrahedralMesh - writing sequential VTK instead");
78  parallelVtk = false;
79  }
80 
81  // If the node ordering flag is set, then we can't do this
82  if (usingOriginalNodeOrdering)
83  {
84  WARNING("Can't write parallel VTK (pvtu) files with original ordering - writing sequential VTK instead");
85  parallelVtk = false;
86  }
87 
88  // Are we now committed to writing .pvtu?
89  if (parallelVtk)
90  {
91  vtk_writer.SetParallelFiles(*pMesh);
92  num_nodes = p_distributed_mesh->GetNumLocalNodes();
93  }
94  }
95 
96  Vec data = p_factory->CreateVec();
97 
98  do // Loop over datasets via MoveOntoNextDataset method in the abstract class
99  {
100  // Make sure that we are never trying to write from an incomplete HDF5 dataset.
101  assert(this->mpReader->GetNumberOfRows() == pMesh->GetNumNodes());
102 
103  unsigned num_timesteps = this->mpReader->GetUnlimitedDimensionValues().size();
104 
105  // Loop over time steps
106  for (unsigned time_step=0; time_step<num_timesteps; time_step++)
107  {
108  // Loop over variables
109  for (unsigned variable=0; variable<this->mNumVariables; variable++)
110  {
111  std::string variable_name = this->mpReader->GetVariableNames()[variable];
112 
113  // Gets variable at this time step from HDF5 archive
114  this->mpReader->GetVariableOverNodes(data, variable_name, time_step);
115 
116  std::vector<double> data_for_vtk;
117  data_for_vtk.resize(num_nodes);
118  std::ostringstream variable_point_data_name;
119  variable_point_data_name << variable_name << "_" << std::setw(6) << std::setfill('0') << time_step;
120 
121  if (parallelVtk)
122  {
123  // Parallel VTU files
124  double *p_data;
125  VecGetArray(data, &p_data);
126  for (unsigned index=0; index<num_nodes; index++)
127  {
128  data_for_vtk[index] = p_data[index];
129  }
130  VecRestoreArray(data, &p_data);
131  }
132  else
133  {
134  // One VTU file
135  ReplicatableVector repl_data(data);
136  for (unsigned index=0; index<num_nodes; index++)
137  {
138  data_for_vtk[index] = repl_data[index];
139  }
140  }
141  // Add this variable into the node "point" data
142  vtk_writer.AddPointData(variable_point_data_name.str(), data_for_vtk);
143  }
144  }
145  }
146  while ( this->MoveOntoNextDataset() );
147 
148  // Tidy up
149  PetscTools::Destroy(data);
150 
151  // Normally the in-memory mesh is converted
152  if (!usingOriginalNodeOrdering)
153  {
154  vtk_writer.WriteFilesUsingMesh(*(this->mpMesh));
155  }
156  else
157  {
158  // In this case we expect the mesh to have been read in from file
160  // Note that the next line will throw if the mesh has not been read from file
161  std::string original_file = this->mpMesh->GetMeshFileBaseName();
162  std::shared_ptr<AbstractMeshReader<ELEMENT_DIM, SPACE_DIM> > p_original_mesh_reader
163  = GenericMeshReader<ELEMENT_DIM, SPACE_DIM>(original_file);
164  vtk_writer.WriteFilesUsingMeshReader(*p_original_mesh_reader);
165  }
166 #endif //CHASTE_VTK
167 }
168 
169 // Explicit instantiation
170 template class Hdf5ToVtkConverter<1,1>;
171 template class Hdf5ToVtkConverter<1,2>;
172 template class Hdf5ToVtkConverter<2,2>;
173 template class Hdf5ToVtkConverter<1,3>;
174 template class Hdf5ToVtkConverter<2,3>;
175 template class Hdf5ToVtkConverter<3,3>;
virtual DistributedVectorFactory * GetDistributedVectorFactory()
Hdf5ToVtkConverter(const FileFinder &rInputDirectory, const std::string &rFileBaseName, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, bool parallelVtk, bool usingOriginalNodeOrdering)
boost::shared_ptr< Hdf5DataReader > mpReader
void SetParallelFiles(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh)
void WriteFilesUsingMeshReader(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)
virtual unsigned GetNumNodes() const
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:352
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
std::string GetRelativePath(const FileFinder &rBasePath) const
Definition: FileFinder.cpp:263
void AddPointData(std::string name, std::vector< double > data)