File: delayed_array.hpp

package info (click to toggle)
r-bioc-alabaster.base 1.6.1%2Bds-2
  • links: PTS, VCS
  • area: main
  • in suites: sid, trixie
  • size: 1,652 kB
  • sloc: cpp: 11,377; sh: 29; makefile: 2
file content (208 lines) | stat: -rw-r--r-- 8,143 bytes parent folder | download | duplicates (2)
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
#ifndef TAKANE_DELAYED_ARRAY_HPP
#define TAKANE_DELAYED_ARRAY_HPP

#include "ritsuko/hdf5/hdf5.hpp"
#include "ritsuko/ritsuko.hpp"
#include "chihaya/chihaya.hpp"

#include "utils_public.hpp"
#include "utils_other.hpp"

#include <vector>
#include <string>
#include <stdexcept>
#include <filesystem>
#include <cstdint>

/**
 * @file delayed_array.hpp
 * @brief Validation for delayed arrays.
 */

namespace takane {

/**
 * @cond
 */
void validate(const std::filesystem::path&, const ObjectMetadata&, Options&);
std::vector<size_t> dimensions(const std::filesystem::path&, const ObjectMetadata&, Options&);
/**
 * @endcond
 */

/**
 * @namespace takane::delayed_array
 * @brief Definitions for delayed arrays.
 */
namespace delayed_array {

/**
 * @cond
 */
namespace internal {

// For efficiency purposes, we just mutate the existing
// 'options.delayed_array_options' rather than making a copy. In this case, we
// need to set 'details_only' to either true or false, depending on whether we
// want to do the full validation; it's important to subsequently reset it back
// to its original setting in the destructor.
struct DetailsOnlyResetter {
    DetailsOnlyResetter(chihaya::Options& o) : options(o), old(options.details_only) {}
    ~DetailsOnlyResetter() {
        options.details_only = old;
    }
private:
    chihaya::Options& options;
    bool old;
};

}
/**
 * @endcond
 */

/**
 * @param path Path to the directory containing a delayed array.
 * @param metadata Metadata for the object, typically read from its `OBJECT` file.
 * @param options Validation options.
 */
inline void validate(const std::filesystem::path& path, const ObjectMetadata& metadata, Options& options) {
    auto vstring = internal_json::extract_version_for_type(metadata.other, "delayed_array");
    auto version = ritsuko::parse_version_string(vstring.c_str(), vstring.size(), /* skip_patch = */ true);
    if (version.major != 1) {
        throw std::runtime_error("unsupported version '" + vstring + "'");
    }

    uint64_t max = 0;
    {
        std::string custom_name = "custom takane seed array";
        auto& custom_options = options.delayed_array_options;
        bool custom_found = (custom_options.array_validate_registry.find(custom_name) != custom_options.array_validate_registry.end());

        // For efficiency purposes, we just mutate the existing
        // 'options.delayed_array_options' rather than making a copy. We need
        // to add a validator for the 'custom takane seed array' type, which
        // checks for valid references to external arrays in 'seeds/'.
        //
        // Note that we respect any existing 'custom takane seed array' setting
        // - possibly from recursive calls to 'delayed_array::validate()', but
        // also possibly from user-provided overrides, in which case we assume
        // that the caller really knows what they're doing.
        //
        // Anyway, all this means that we only mutate chihaya::Options if there
        // is no existing custom takane function. However, if we do so, we need
        // to restore the original state before function exit, hence the
        // destructor for RAII-based clean-up.
        struct ValidateResetter {
            ValidateResetter(chihaya::Options& o, const std::string& n, bool f) : options(o), name(n), found(f) {}
            ~ValidateResetter() {
                if (!found) {
                    options.array_validate_registry.erase(name);
                }
            }
        private:
            chihaya::Options& options;
            const std::string& name;
            bool found; 
        };
        [[maybe_unused]] ValidateResetter v(custom_options, custom_name, custom_found);

        if (!custom_found) {
            custom_options.array_validate_registry[custom_name] = [&](const H5::Group& handle, const ritsuko::Version& version, chihaya::Options& ch_options) -> chihaya::ArrayDetails {
                auto details = chihaya::custom_array::validate(handle, version, ch_options);

                auto dhandle = ritsuko::hdf5::open_dataset(handle, "index");
                if (ritsuko::hdf5::exceeds_integer_limit(dhandle, 64, false)) {
                    throw std::runtime_error("'index' should have a datatype that fits into a 64-bit unsigned integer");
                }

                auto index = ritsuko::hdf5::load_scalar_numeric_dataset<uint64_t>(dhandle);
                auto seed_path = path / "seeds" / std::to_string(index);
                auto seed_meta = read_object_metadata(seed_path);
                ::takane::validate(seed_path, seed_meta, options);

                auto seed_dims = ::takane::dimensions(seed_path, seed_meta, options);
                if (seed_dims.size() != details.dimensions.size()) {
                    throw std::runtime_error("dimensionality of 'seeds/" + std::to_string(index) + "' is not consistent with 'dimensions'");
                }

                for (size_t d = 0, ndims = seed_dims.size(); d < ndims; ++d) {
                    if (seed_dims[d] != details.dimensions[d]) {
                        throw std::runtime_error("dimension extents of 'seeds/" + std::to_string(index) + "' is not consistent with 'dimensions'");
                    }
                }

                if (index >= max) {
                    max = index + 1;
                }
                return details;
            };
        }

        auto apath = path / "array.h5";
        auto fhandle = ritsuko::hdf5::open_file(apath);
        auto ghandle = ritsuko::hdf5::open_group(fhandle, "delayed_array");
        ritsuko::Version chihaya_version = chihaya::extract_version(ghandle);
        if (chihaya_version.lt(1, 1, 0)) {
            throw std::runtime_error("version of the chihaya specification should be no less than 1.1");
        }

        // Again, using RAII to reset the 'details_only' flag to its original
        // state after we're done with it.
        [[maybe_unused]] internal::DetailsOnlyResetter o(custom_options);
        custom_options.details_only = false;

        chihaya::validate(ghandle, chihaya_version, custom_options);
    }

    size_t found = 0;
    auto seed_path = path / "seeds";
    if (std::filesystem::exists(seed_path)) {
        found = internal_other::count_directory_entries(seed_path);
    }
    if (max != found) {
        throw std::runtime_error("number of objects in 'seeds' is not consistent with the number of 'index' references in 'array.h5'");
    }
}

/**
 * @param path Path to the directory containing a delayed array.
 * @param metadata Metadata for the object, typically read from its `OBJECT` file.
 * @param options Validation options.
 * @return Extent of the first dimension.
 */
inline size_t height(const std::filesystem::path& path, [[maybe_unused]] const ObjectMetadata& metadata, Options& options) {
    auto& chihaya_options = options.delayed_array_options;
    [[maybe_unused]] internal::DetailsOnlyResetter o(chihaya_options);
    chihaya_options.details_only = true;

    auto apath = path / "array.h5";
    auto fhandle = ritsuko::hdf5::open_file(apath);
    auto ghandle = ritsuko::hdf5::open_group(fhandle, "delayed_array");
    auto output = chihaya::validate(ghandle, chihaya_options);
    return output.dimensions[0];
}

/**
 * @param path Path to the directory containing a delayed array.
 * @param metadata Metadata for the object, typically read from its `OBJECT` file.
 * @param options Validation options.
 * @return Dimensions of the array.
 */
inline std::vector<size_t> dimensions(const std::filesystem::path& path, [[maybe_unused]] const ObjectMetadata& metadata, Options& options) {
    auto& chihaya_options = options.delayed_array_options;
    [[maybe_unused]] internal::DetailsOnlyResetter o(chihaya_options);
    chihaya_options.details_only = true;

    auto apath = path / "array.h5";
    auto fhandle = ritsuko::hdf5::open_file(apath);
    auto ghandle = ritsuko::hdf5::open_group(fhandle, "delayed_array");
    auto output = chihaya::validate(ghandle, chihaya_options);
    return std::vector<size_t>(output.dimensions.begin(), output.dimensions.end());
}

}

}

#endif