Coverage Report

Created: 2025-07-09 06:49

/src/htslib/test/fuzz/hts_open_fuzzer.c
Line
Count
Source (jump to first uncovered line)
1
/*  test/fuzz/hts_open_fuzzer.c -- Fuzz driver for hts_open.
2
3
    Copyright (C) 2018 Google LLC.
4
    Copyright (C) 2019-2020, 2023 Genome Research Ltd.
5
6
    Author: Markus Kusano <kusano@google.com>
7
8
Permission is hereby granted, free of charge, to any person obtaining a copy
9
of this software and associated documentation files (the "Software"), to deal
10
in the Software without restriction, including without limitation the rights
11
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
12
copies of the Software, and to permit persons to whom the Software is
13
furnished to do so, subject to the following conditions:
14
15
The above copyright notice and this permission notice shall be included in
16
all copies or substantial portions of the Software.
17
18
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
19
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
20
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
21
THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
22
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
23
FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
24
DEALINGS IN THE SOFTWARE.  */
25
26
#include <config.h>
27
28
#include <stddef.h>
29
#include <stdint.h>
30
#include <stdio.h>
31
#include <stdlib.h>
32
#include <string.h>
33
#include <unistd.h>
34
35
#include "../../htslib/hfile.h"
36
#include "../../htslib/hts.h"
37
#include "../../htslib/sam.h"
38
#include "../../htslib/vcf.h"
39
40
23.5k
static void hts_close_or_abort(htsFile* file) {
41
23.5k
    if (hts_close(file) != 0) {
42
0
        abort();
43
0
    }
44
23.5k
}
45
46
static void view_sam(const uint8_t *data, size_t size, char *mode,
47
28.1k
                     int close_abort) {
48
28.1k
    uint8_t *copy = malloc(size);
49
28.1k
    if (copy == NULL) {
50
0
        abort();
51
0
    }
52
28.1k
    memcpy(copy, data, size);
53
54
28.1k
    hFILE *memfile = hopen("mem:", "rb:", copy, size);
55
28.1k
    if (memfile == NULL) {
56
0
        free(copy);
57
0
        return;
58
0
    }
59
60
28.1k
    htsFile *in = hts_hopen(memfile, "data", "rb");
61
28.1k
    if (in == NULL) {
62
0
        if (hclose(memfile) != 0)
63
0
            abort();
64
0
        return;
65
0
    }
66
67
28.1k
    samFile *out = sam_open("/dev/null", mode);
68
28.1k
    if (!out)
69
0
        abort();
70
71
#ifdef FUZZ_FAI
72
    // Not critical if this doesn't work, but can test more if
73
    // we're in the right location.
74
    //
75
    // We can't rely on what the pwd is for the OSS-fuzz so we don't enable
76
    // this by default.
77
    if (hts_set_fai_filename(out, "../c2.fa") < 0) {
78
        static int warned = 0;
79
        if (!warned) {
80
            warned = 1;
81
            fprintf(stderr, "Warning couldn't find the c2.fa file\n");
82
        }
83
    }
84
#endif
85
86
28.1k
    sam_hdr_t *hdr = sam_hdr_read(in);
87
28.1k
    if (hdr == NULL) {
88
1.17k
        if (close_abort)
89
786
            hts_close_or_abort(out);
90
393
        else
91
393
            hts_close(out);
92
1.17k
        hts_close(in);
93
1.17k
        return;
94
1.17k
    }
95
96
    // This will force the header to be parsed.
97
26.9k
    (void) sam_hdr_count_lines(hdr, "SQ");
98
99
26.9k
    if (sam_hdr_write(out, hdr) != 0)
100
1.05k
        goto err;
101
102
25.9k
    bam1_t *b = bam_init1();
103
25.9k
    if (b == NULL)
104
0
        goto err;
105
106
4.58M
    while (sam_read1(in, hdr, b) >= 0) {
107
4.55M
        if (sam_write1(out, hdr, b) < 0)
108
804
            break;
109
4.55M
    }
110
25.9k
    bam_destroy1(b);
111
112
26.9k
 err:
113
26.9k
    sam_hdr_destroy(hdr);
114
26.9k
    if (close_abort)
115
17.9k
        hts_close_or_abort(out);
116
8.98k
    else
117
8.98k
        hts_close(out);
118
26.9k
    hts_close(in);
119
26.9k
}
120
121
4.75k
static void view_vcf(const uint8_t *data, size_t size, char *mode) {
122
4.75k
    uint8_t *copy = malloc(size);
123
4.75k
    if (copy == NULL) {
124
0
        abort();
125
0
    }
126
4.75k
    memcpy(copy, data, size);
127
128
4.75k
    hFILE *memfile = hopen("mem:", "rb:", copy, size);
129
4.75k
    if (memfile == NULL) {
130
0
        free(copy);
131
0
        return;
132
0
    }
133
134
4.75k
    htsFile *in = hts_hopen(memfile, "data", "rb");
135
4.75k
    if (in == NULL) {
136
0
        if (hclose(memfile) != 0)
137
0
            abort();
138
0
        return;
139
0
    }
140
141
4.75k
    vcfFile *out = vcf_open("/dev/null", mode);
142
4.75k
    if (!out)
143
0
        abort();
144
145
4.75k
    bcf_hdr_t *hdr = bcf_hdr_read(in);
146
4.75k
    if (hdr == NULL) {
147
1.71k
        hts_close_or_abort(out);
148
1.71k
        hts_close(in);
149
1.71k
        return;
150
1.71k
    }
151
152
3.03k
    if (bcf_hdr_write(out, hdr) != 0)
153
0
        goto err;
154
155
3.03k
    bcf1_t *rec = bcf_init();
156
3.03k
    if (rec == NULL)
157
0
        goto err;
158
159
19.9k
    while (bcf_read(in, hdr, rec) >= 0) {
160
18.1k
        if (bcf_write(out, hdr, rec) < 0)
161
1.24k
            break;
162
18.1k
    }
163
3.03k
    bcf_destroy(rec);
164
165
3.03k
 err:
166
3.03k
    bcf_hdr_destroy(hdr);
167
3.03k
    hts_close_or_abort(out);
168
3.03k
    hts_close(in);
169
3.03k
}
170
171
14.2k
int LLVMFuzzerTestOneInput(const uint8_t *data, size_t size) {
172
    // Only data as a mem file purely for purposes of determining format
173
14.2k
    hFILE *memfile;
174
14.2k
    uint8_t *copy = malloc(size);
175
14.2k
    if (copy == NULL) {
176
0
        abort();
177
0
    }
178
14.2k
    memcpy(copy, data, size);
179
    // hopen does not take ownership of `copy`, but hts_hopen does.
180
14.2k
    memfile = hopen("mem:", "rb:", copy, size);
181
14.2k
    if (memfile == NULL) {
182
0
        free(copy);
183
0
        return 0;
184
0
    }
185
186
14.2k
    htsFile *ht_file = hts_hopen(memfile, "data", "rb");
187
14.2k
    if (ht_file == NULL) {
188
2.24k
        if (hclose(memfile) != 0) {
189
0
            abort();
190
0
        }
191
2.24k
        return 0;
192
2.24k
    }
193
11.9k
    int ftype = ht_file->format.category;
194
11.9k
    hts_close(ht_file);
195
196
    // Now repeat a read-write loop multiple times per input, testing
197
    // encoding in all output formats.
198
    // (Although we could just ignore ftype and do all 5 for all inputs)
199
11.9k
    switch (ftype) {
200
9.37k
        case sequence_data:
201
9.37k
            view_sam(data, size, "w",  1); // SAM
202
9.37k
            view_sam(data, size, "wb", 1); // BAM
203
9.37k
            view_sam(data, size, "wc", 0); // CRAM
204
9.37k
            break;
205
2.37k
        case variant_data:
206
2.37k
            view_vcf(data, size, "w");     // VCF
207
2.37k
            view_vcf(data, size, "wb");    // BCF
208
2.37k
            break;
209
244
        default:
210
244
            break;
211
11.9k
    }
212
11.9k
    return 0;
213
11.9k
}