File: grib_get_keys.c

package info (click to toggle)
eccodes 2.20.0-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 400,332 kB
  • sloc: ansic: 167,977; makefile: 21,348; sh: 10,719; f90: 5,927; python: 4,831; perl: 3,031; javascript: 1,427; yacc: 818; lex: 356; awk: 66
file content (147 lines) | stat: -rw-r--r-- 5,130 bytes parent folder | download
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
/*
 * (C) Copyright 2005- ECMWF.
 *
 * This software is licensed under the terms of the Apache Licence Version 2.0
 * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
 *
 * In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
 * virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
 */

/*
 * C Implementation: grib_get_keys
 *
 * Description: how to get values using keys from GRIB messages
 *
 */
#include <stdio.h>
#include <stdlib.h>
#include <assert.h>
#include "eccodes.h"

int main(int argc, char** argv)
{
    int err           = 0;
    double* values    = NULL;
    size_t values_len = 0;
    size_t i = 0, len = 0;

    double latitudeOfFirstGridPointInDegrees;
    double longitudeOfFirstGridPointInDegrees;
    double latitudeOfLastGridPointInDegrees;
    double longitudeOfLastGridPointInDegrees;

    double jDirectionIncrementInDegrees;
    double iDirectionIncrementInDegrees;

    long numberOfPointsAlongAParallel;
    long numberOfPointsAlongAMeridian;

    double average    = 0;
    char* packingType = NULL;

    FILE* in             = NULL;
    const char* filename = "../../data/regular_latlon_surface.grib1";
    codes_handle* h      = NULL;

    in = fopen(filename, "rb");
    if (!in) {
        fprintf(stderr, "Error: unable to open file %s\n", filename);
        return 1;
    }

    /* create new handle from the first message in the file*/
    h = codes_handle_new_from_file(0, in, PRODUCT_GRIB, &err);
    if (h == NULL) {
        fprintf(stderr, "Error: unable to create handle from file %s\n", filename);
        return 1;
    }
    fclose(in);

    /* Store the filename in the key "file" for this handle */
    len = strlen(filename);
    CODES_CHECK(codes_set_string(h, "file", filename, &len), 0);

    /* get as a long*/
    CODES_CHECK(codes_get_long(h, "Ni", &numberOfPointsAlongAParallel), 0);
    printf("numberOfPointsAlongAParallel=%ld\n", numberOfPointsAlongAParallel);

    /* get as a long*/
    CODES_CHECK(codes_get_long(h, "Nj", &numberOfPointsAlongAMeridian), 0);
    printf("numberOfPointsAlongAMeridian=%ld\n", numberOfPointsAlongAMeridian);

    /* get as a double*/
    CODES_CHECK(codes_get_double(h, "latitudeOfFirstGridPointInDegrees", &latitudeOfFirstGridPointInDegrees), 0);
    printf("latitudeOfFirstGridPointInDegrees=%g\n", latitudeOfFirstGridPointInDegrees);

    /* get as a double*/
    CODES_CHECK(codes_get_double(h, "longitudeOfFirstGridPointInDegrees", &longitudeOfFirstGridPointInDegrees), 0);
    printf("longitudeOfFirstGridPointInDegrees=%g\n", longitudeOfFirstGridPointInDegrees);

    /* get as a double*/
    CODES_CHECK(codes_get_double(h, "latitudeOfLastGridPointInDegrees", &latitudeOfLastGridPointInDegrees), 0);
    printf("latitudeOfLastGridPointInDegrees=%g\n", latitudeOfLastGridPointInDegrees);

    /* get as a double*/
    CODES_CHECK(codes_get_double(h, "longitudeOfLastGridPointInDegrees", &longitudeOfLastGridPointInDegrees), 0);
    printf("longitudeOfLastGridPointInDegrees=%g\n", longitudeOfLastGridPointInDegrees);

    /* get as a double*/
    CODES_CHECK(codes_get_double(h, "jDirectionIncrementInDegrees", &jDirectionIncrementInDegrees), 0);
    printf("jDirectionIncrementInDegrees=%g\n", jDirectionIncrementInDegrees);

    /* get as a double*/
    CODES_CHECK(codes_get_double(h, "iDirectionIncrementInDegrees", &iDirectionIncrementInDegrees), 0);
    printf("iDirectionIncrementInDegrees=%g\n", iDirectionIncrementInDegrees);

    /* get as string */
    CODES_CHECK(codes_get_length(h, "packingType", &len), 0);
    packingType = (char*)malloc(len * sizeof(char));
    codes_get_string(h, "packingType", packingType, &len);
    printf("packingType=%s\n", packingType);
    free(packingType);

    /* get the size of the values array*/
    CODES_CHECK(codes_get_size(h, "values", &values_len), 0);

    values = (double*)malloc(values_len * sizeof(double));

    /* get data values*/
    CODES_CHECK(codes_get_double_array(h, "values", values, &values_len), 0);

    average = 0;
    for (i = 0; i < values_len; i++)
        average += values[i];

    average /= (double)values_len;

    free(values);

    printf("There are %d values, average is %g\n", (int)values_len, average);

    {
        int eq = 0;
        /* Now retrieve the value of the key "file" */
        char file[256] = {0,};
        CODES_CHECK(codes_get_length(h, "file", &len), 0);
        assert(len == 1 + strlen(filename));
        codes_get_string(h, "file", file, &len);
        eq = strcmp(file, filename);
        if (eq != 0) assert(!"file and filename not equal");
    }

    {
        /* Example of getting bytes */
        const char* name        = "reservedNeedNotBePresent";
        unsigned char* byte_val = NULL;
        size_t keySize          = 0;
        CODES_CHECK(codes_get_size(h, name, &keySize), 0);
        byte_val = (unsigned char*)malloc(keySize * sizeof(unsigned char));
        GRIB_CHECK(codes_get_bytes(h, name, byte_val, &keySize), name);
        free(byte_val);
    }

    codes_handle_delete(h);

    return 0;
}