This file is indexed.

/usr/include/dolfin/io/HDF5File.h is in libdolfin-dev 2016.2.0-2.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
// Copyright (C) 2012 Chris N. Richardson
//
// This file is part of DOLFIN.
//
// DOLFIN is free software: you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// DOLFIN is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
//
// Modified by Garth N. Wells, 2012

#ifndef __DOLFIN_HDF5FILE_H
#define __DOLFIN_HDF5FILE_H

#ifdef HAS_HDF5

#include <string>
#include <utility>
#include <vector>

#include <dolfin/common/MPI.h>
#include <dolfin/common/Variable.h>
#include <dolfin/geometry/Point.h>
#include "HDF5Attribute.h"
#include "HDF5Interface.h"

namespace dolfin
{

  class CellType;
  class Function;
  class GenericVector;
  class LocalMeshData;
  class Mesh;
  template<typename T> class MeshFunction;
  template<typename T> class MeshValueCollection;
  class HDF5Attribute;

  class HDF5File : public Variable
  {

  public:

    /// Constructor. file_mode should be "a" (append),
    /// "w" (write) or "r" (read).
    HDF5File(MPI_Comm comm, const std::string filename,
             const std::string file_mode);

    /// Destructor
    ~HDF5File();

    /// Close file
    void close();

    /// Flush buffered I/O to disk
    void flush();

    /// Write points to file
    void write(const std::vector<Point>& points, const std::string name);

    /// Write simple vector of double to file
    void write(const std::vector<double>& values, const std::string name);

    /// Write Vector to file in a format suitable for re-reading
    void write(const GenericVector& x, const std::string name);

    /// Read vector from file and optionally re-use any partitioning
    /// that is available in the file
    void read(GenericVector& x, const std::string dataset_name,
              const bool use_partition_from_file) const;

    /// Write Mesh to file in a format suitable for re-reading
    void write(const Mesh& mesh, const std::string name);

    /// Write Mesh of given cell dimension to file in a format
    /// suitable for re-reading
    void write(const Mesh& mesh, const std::size_t cell_dim,
               const std::string name);

    /// Write Function to file in a format suitable for re-reading
    void write(const Function& u, const std::string name);

    /// Write Function to file with a timestamp
    void write(const Function& u, const std::string name, double timestamp);

    /// Read Function from file and distribute data according to the
    /// Mesh and dofmap associated with the Function.  If the 'name'
    /// refers to a HDF5 group, then it is assumed that the Function
    /// data is stored in the datasets within that group.  If the
    /// 'name' refers to a HDF5 dataset within a group, then it is
    /// assumed that it is a Vector, and the Function will be filled
    /// from that Vector
    void read(Function& u, const std::string name);

    /// Read Mesh from file, using attribute data (e.g., cell type) stored
    /// in the HDF5 file. Optionally re-use any partition data
    /// in the file. This function requires all necessary data for
    /// constructing a Mesh to be present in the HDF5 file.
    void read(Mesh& mesh, const std::string data_path,
              bool use_partition_from_file) const;

    /// Construct Mesh with paths to topology and geometry datasets, and
    /// providing essential meta-data, e.g. geometric dimension and
    /// cell type. If this data is available in the HDF5 file, it will
    /// be checked for consistency. Set expected_num_global_cells to a
    /// negative value if not known.
    ///
    /// This function is typically called when using the XDMF format,
    /// in which case the meta data has alreayd been read from an XML
    /// file
    void read(Mesh& input_mesh,
              const std::string topology_path,
              const std::string geometry_path,
              const int gdim , const CellType& cell_type,
              const std::int64_t expected_num_global_cells,
              const std::int64_t expected_num_global_points,
              bool use_partition_from_file) const;

    /// Write MeshFunction to file in a format suitable for re-reading
    void write(const MeshFunction<std::size_t>& meshfunction,
               const std::string name);

    /// Write MeshFunction to file in a format suitable for re-reading
    void write(const MeshFunction<int>& meshfunction, const std::string name);

    /// Write MeshFunction to file in a format suitable for re-reading
    void write(const MeshFunction<double>& meshfunction,
               const std::string name);

    /// Write MeshFunction to file in a format suitable for re-reading
    void write(const MeshFunction<bool>& meshfunction, const std::string name);

    /// Read MeshFunction from file
    void read(MeshFunction<std::size_t>& meshfunction,
              const std::string name) const;

    /// Read MeshFunction from file
    void read(MeshFunction<int>& meshfunction, const std::string name) const;

    /// Read MeshFunction from file
    void read(MeshFunction<double>& meshfunction,
              const std::string name) const;

    /// Read MeshFunction from file
    void read(MeshFunction<bool>& meshfunction, const std::string name) const;

    /// Write MeshValueCollection to file
    void write(const MeshValueCollection<std::size_t>& mesh_values,
               const std::string name);

    /// Write MeshValueCollection to file
    void write(const MeshValueCollection<double>& mesh_values,
               const std::string name);

    /// Write MeshValueCollection to file
    void write(const MeshValueCollection<bool>& mesh_values,
               const std::string name);

    /// Read MeshValueCollection from file
    void read(MeshValueCollection<std::size_t>& mesh_values,
              const std::string name) const;

    /// Read MeshValueCollection from file
    void read(MeshValueCollection<double>& mesh_values,
              const std::string name) const;

    /// Read MeshValueCollection from file
    void read(MeshValueCollection<bool>& mesh_values,
              const std::string name) const;

    /// Check if dataset exists in HDF5 file
    bool has_dataset(const std::string dataset_name) const;

    // Get/set attributes of an existing dataset
    HDF5Attribute attributes(const std::string dataset_name);

    /// Set the MPI atomicity
    void set_mpi_atomicity(bool atomic);

    /// Get the MPI atomicity
    bool get_mpi_atomicity() const;

    hid_t h5_id() const
    { return _hdf5_file_id; }

  private:

    // Friend
    friend class XDMFFile;
    friend class TimeSeries;

    // Write a MeshFunction to file
    template <typename T>
    void write_mesh_function(const MeshFunction<T>& meshfunction,
                             const std::string name);

    // Read a MeshFunction from file
    template <typename T>
    void read_mesh_function(MeshFunction<T>& meshfunction,
                            const std::string name) const;

    // Write a MeshValueCollection to file (old format)
    template <typename T>
    void write_mesh_value_collection_old(
                                     const MeshValueCollection<T>& mesh_values,
                                     const std::string name);

    // Write a MeshValueCollection to file (new version using vertex indices)
    template <typename T>
    void write_mesh_value_collection(const MeshValueCollection<T>& mesh_values,
                                     const std::string name);

    // Read a MeshValueCollection from file
    template <typename T>
    void read_mesh_value_collection(MeshValueCollection<T>& mesh_values,
                                    const std::string name) const;

    // Read a MeshValueCollection (old format)
    template <typename T>
    void read_mesh_value_collection_old(MeshValueCollection<T>& mesh_values,
                                    const std::string name) const;

    // Write contiguous data to HDF5 data set. Data is flattened into
    // a 1D array, e.g. [x0, y0, z0, x1, y1, z1] for a vector in 3D
    template <typename T>
    void write_data(const std::string dataset_name,
                    const std::vector<T>& data,
                    const std::vector<std::int64_t> global_size,
                    bool use_mpi_io);

    // HDF5 file descriptor/handle
    hid_t _hdf5_file_id;

    // MPI communicator
    MPI_Comm _mpi_comm;
  };

  //---------------------------------------------------------------------------
  // Needs to go here, because of use in XDMFFile.cpp
  template <typename T>
  void HDF5File::write_data(const std::string dataset_name,
                            const std::vector<T>& data,
                            const std::vector<std::int64_t> global_size,
                            bool use_mpi_io)
  {
    dolfin_assert(_hdf5_file_id > 0);
    dolfin_assert(global_size.size() > 0);

    // Get number of 'items'
    std::size_t num_local_items = 1;
    for (std::size_t i = 1; i < global_size.size(); ++i)
      num_local_items *= global_size[i];
    num_local_items = data.size()/num_local_items;

    // Compute offset
    const std::size_t offset = MPI::global_offset(_mpi_comm, num_local_items,
                                                  true);
    std::pair<std::size_t, std::size_t> range(offset,
                                              offset + num_local_items);

    // Write data to HDF5 file
    const bool chunking = parameters["chunking"];
    // Ensure dataset starts with '/'
    std::string dset_name(dataset_name);
    if (dset_name[0] != '/')
      dset_name = "/" + dataset_name;

    HDF5Interface::write_dataset(_hdf5_file_id, dset_name, data,
                                 range, global_size, use_mpi_io, chunking);
  }
  //---------------------------------------------------------------------------

}

#endif
#endif