File: bufr_read_scatterometer.c

package info (click to toggle)
eccodes 2.44.2-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 150,256 kB
  • sloc: cpp: 163,056; ansic: 26,308; sh: 21,602; f90: 6,854; perl: 6,363; python: 5,087; java: 2,226; javascript: 1,427; yacc: 854; fortran: 543; lex: 359; makefile: 274; xml: 183; awk: 66
file content (138 lines) | stat: -rw-r--r-- 4,785 bytes parent folder | download | duplicates (3)
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
/*
 * (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: bufr_read_scatterometer
 *
 * Description: how to read data for a given beam from scatterometer BUFR messages.
 *
 */

/*
 * Please note that scatterometer data can be encoded in various ways in BUFR. Therefore the code
 * below might not work directly for other types of messages than the one used in the
 * example. It is advised to use bufr_dump to understand the structure of the messages.
 */

#include "eccodes.h"

int main(int argc, char* argv[])
{
    FILE* in = NULL;

    /* Message handle. Required in all the eccodes calls acting on a message.*/
    codes_handle* h = NULL;

    double *lat = NULL, *lon = NULL, *bscatter = NULL;
    long numObs = 0;
    size_t len  = 0;
    int i, err = 0;
    int cnt            = 0;
    const char* infile = "../../data/bufr/asca_139.bufr";
    char key_name[128];

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

    /* Loop over the messages in the BUFR file */
    while ((h = codes_handle_new_from_file(NULL, in, PRODUCT_BUFR, &err)) != NULL || err != CODES_SUCCESS) {
        if (h == NULL) {
            fprintf(stderr, "Error: unable to create handle for message %d\n", cnt);
            cnt++;
            continue;
        }

        printf("message: %d\n", cnt);

        /* We need to instruct ecCodes to expand the descriptors
         * i.e. unpack the data values */
        CODES_CHECK(codes_set_long(h, "unpack", 1), 0);

        /* The BUFR file contains a single message with 2016 subsets in a compressed form.
         * It means each subset has exactly the same structure: they store one location with
         * several beams and one backscatter value in each beam.
         *
         * To print the backScatter values for beamIdentifier=2 from all the subsets
         * we will simply access the key by condition (see below) */

        /* Get the total number of subsets. */
        CODES_CHECK(codes_get_long(h, "numberOfSubsets", &numObs), 0);

        printf("Number of values: %ld\n", numObs);

        /* Get latitude */
        snprintf(key_name, sizeof(key_name), "latitude");

        /* Check the size (including all the subsets) */
        CODES_CHECK(codes_get_size(h, key_name, &len), 0);
        if (len != numObs) {
            fprintf(stderr, "Error: inconsistent number of %s values found!\n", key_name);
            return 1;
        }

        /* Allocate memory for the values to be read. Each
         * parameter must have the same number of values. */
        lat = (double*)malloc(numObs * sizeof(double));

        /* Get the values (from all the subsets) */
        CODES_CHECK(codes_get_double_array(h, key_name, lat, &len), 0);

        /* Get longitude */
        snprintf(key_name, sizeof(key_name), "longitude");

        /* Check the size (including all the subsets) */
        CODES_CHECK(codes_get_size(h, key_name, &len), 0);
        if (len != numObs) {
            fprintf(stderr, "Error: inconsistent number of %s values found!\n", key_name);
            return 1;
        }

        /* Get the values (from all the subsets) */
        lon = (double*)malloc(numObs * sizeof(double));
        CODES_CHECK(codes_get_double_array(h, key_name, lon, &len), 0);

        /* Get backScatter for beam two. We use an access by condition for this key. */
        snprintf(key_name, sizeof(key_name), "/beamIdentifier=2/backscatter");

        /* Check the size (including all the subsets) */
        CODES_CHECK(codes_get_size(h, key_name, &len), 0);
        if (len != numObs) {
            fprintf(stderr, "Error: inconsistent number of %s values found!\n", key_name);
            return 1;
        }

        /* Get the values (from all the subsets) */
        bscatter = (double*)malloc(numObs * sizeof(double));
        CODES_CHECK(codes_get_double_array(h, key_name, bscatter, &len), 0);

        /* Print the values */
        printf("pixel   lat    lon     backscatter    \n");
        printf("-------------------------------\n");
        for (i = 0; i < numObs; i++) {
            printf("%4d %.3f %.3f %.3f \n", i + 1, lat[i], lon[i], bscatter[i]);
        }

        /* Delete handle */
        codes_handle_delete(h);

        /* Release memory */
        free(lat);
        free(lon);
        free(bscatter);

        cnt++;
    }

    fclose(in);
    return 0;
}