IODA
05-ObsGroup.py
Go to the documentation of this file.
1 #
2 # (C) Copyright 2020-2021 UCAR
3 #
4 # This software is licensed under the terms of the Apache Licence Version 2.0
5 # which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6 #
7 
8 # The ObsGroup class is derived from the Group class and provides some help in
9 # organizing your groups, variables, attributes and dimension scales into a cohesive
10 # structure intended to house observation data. In this case "structure" refers to the
11 # hierarchical layout of the groups and the proper management of dimension scales
12 # associated with the variables.
13 #
14 # The ObsGroup and underlying layout policies (internal to ioda-engines) present a stable
15 # logical group hierarchical layout to the client while keeping the actual layout implemented
16 # in the backend open to change. The logical "frontend" layout appears to the client to be
17 # as shown below:
18 #
19 # layout notes
20 #
21 # / top-level group
22 # nlocs dimension scales (variables, coordinate values)
23 # nchans
24 # ...
25 # ObsValue/ group: observational measurement values
26 # brightness_temperature variable: Tb, 2D, nlocs X nchans
27 # air_temperature variable: T, 1D, nlocs
28 # ...
29 # ObsError/ group: observational error estimates
30 # brightness_temperature
31 # air_temperature
32 # ...
33 # PreQC/ group: observational QC marks from data provider
34 # brightness_temperature
35 # air_temperature
36 # ...
37 # MetaData/ group: meta data associated with locations
38 # latitude
39 # longitude
40 # datetime
41 # ...
42 # ...
43 #
44 # We intend to keep this layout stable so that the client interface remains stable.
45 # The actual layout used in the various backends can optionally be organized differently
46 # according to their needs.
47 #
48 # The ObsGroup class also assists with the management of dimension scales. For example, if
49 # a dimension is resized, the ObsGroup::resize function will resize the dimension scale
50 # along with all variables that use that dimension scale.
51 #
52 # The basic ideas is to dimension observation data with nlocs as the first dimension, and
53 # allow nlocs to be resizable so that it's possible to incrementally append data along
54 # the nlocs (1st) dimension. For data that have rank > 1, the second through nth dimensions
55 # are of fixed size. For example, brightness_temperature can be store as 2D data with
56 # dimensions (nlocs, nchans).
57 
58 import os
59 import sys
60 
61 import ioda
62 import numpy as np
63 
65  name = "Example-05a-python.hdf5",
66  mode = ioda.Engines.BackendCreateModes.Truncate_If_Exists)
67 
68 # We have opened the file, but now we want to turn it into an ObsGroup.
69 # To do this, we need to first specify our dimensions. Then, we call
70 # the ObsGroup.generate function.
71 
72 
73 numLocs = 40
74 numChans = 30
75 
76 # This is a list of the dimensions that we want in our ObsGroup.
77 # Both dimensions can be represented as 32-bit integers.
78 # They are called 'nlocs' and 'nchans'.
79 #
80 # The NewDimensionScale.* functions take four parameters:
81 # 1. The name of the dimension scale.
82 # 2. The initial length (size) of this scale. If you have 40 locations,
83 # then the nlocs scale's initial size should be 40.
84 # 3. The maximum length of the scale. By setting this to ioda.Unlimited,
85 # you can resize all of the variables depending on a scale and
86 # append new data without recreating the ObsSpace. This doesn't make
87 # sense for all possible dimensions... instrument channels are usually
88 # fixed, so the max length should match the initial length in those cases.
89 # 4. A "hint" for the chunking size used by Variables that depend on this
90 # scale. When reading data from a disk or in memory, ioda stores data in
91 # "blocks". The chunking size determines the size of these blocks.
92 # For details on block sizes see https://support.hdfgroup.org/HDF5/doc/Advanced/Chunking/
93 #
94 # Let's make a list of the new dimension scales.
95 newDims = [ioda.NewDimensionScale.int32('nlocs', numLocs, ioda.Unlimited, numLocs),
96  ioda.NewDimensionScale.int32('nchans', numChans, numChans, numChans)]
97 
98 # ObsGroup.generate takes a Group argument (the backend we just created
99 # above) and a list of dimension creation scales. It then creates the dimensions
100 # and performs initialization of an ObsGroup.
101 og = ioda.ObsGroup.generate(g, newDims)
102 
103 # You can open the scales using the names that you provided in the
104 # NewDimensionScale calls.
105 nlocsVar = og.vars.open('nlocs')
106 nchansVar = og.vars.open('nchans')
107 
108 # Just setting some sensible defaults: compress the stored data and use
109 # a fill value of -999.
111 p1.compressWithGZIP()
112 p1.setFillValue.float(-999)
113 
114 # Next let's create the variables. The variable names should be specified using the
115 # hierarchy as described above. For example, the variable brightness_temperature
116 # in the group ObsValue is specified in a string as "ObsValue/brightness_temperature".
117 tbName = "ObsValue/brightness_temperature"
118 latName = "MetaData/latitude"
119 lonName = "MetaData/longitude"
120 
121 # We create three variables to store brightness temperature, latitude and longitude.
122 # We attach the appropriate scales with the scales option.
123 tbVar = g.vars.create(tbName, ioda.Types.float, scales=[nlocsVar, nchansVar], params=p1)
124 latVar = g.vars.create(latName, ioda.Types.float, scales=[nlocsVar], params=p1)
125 lonVar = g.vars.create(lonName, ioda.Types.float, scales=[nlocsVar], params=p1)
126 
127 # Let's set some attributes on the variables to describe what they mean,
128 # their ranges, their units, et cetera.
129 tbVar.atts.create("coordinates", ioda.Types.str, [1]).writeDatum.str("longitude latitude nchans")
130 tbVar.atts.create("long_name", ioda.Types.str, [1]).writeDatum.str("fictional brightness temperature")
131 tbVar.atts.create("units", ioda.Types.str, [1]).writeDatum.str("K")
132 tbVar.atts.create("valid_range", ioda.Types.float, [2]).writeVector.float([100.0, 400.0])
133 
134 latVar.atts.create("long_name", ioda.Types.str, [1]).writeDatum.str("latitude")
135 latVar.atts.create("units", ioda.Types.str, [1]).writeDatum.str("degrees_north")
136 latVar.atts.create("valid_range", ioda.Types.float, [2]).writeVector.float([-90, 90])
137 
138 lonVar.atts.create("long_name", ioda.Types.str, [1]).writeDatum.str("degrees_east")
139 lonVar.atts.create("units", ioda.Types.str, [1]).writeDatum.str("degrees_north")
140 lonVar.atts.create("valid_range", ioda.Types.float, [2]).writeVector.float([-360, 360])
141 
142 # Let's create some data for this example.
143 # We are generating filler data.
144 # Perhaps consult numpy.fromfunction's documentation.
145 
146 midLoc = numLocs / 2
147 midChan = numChans / 2
148 
149 # Latitude and longitude are 1-D numpy arrays, with length numLocs.
150 # Each array index is fed into a lambda function that generates a value
151 # at each location.
152 lonData = np.fromfunction(lambda i: 3 * (i % 8), (numLocs,), dtype=float)
153 latData = np.fromfunction(lambda i: 3 * (i // 8), (numLocs,), dtype=float)
154 # Brightness temperatures are two-dimensional. They have dimensions of
155 # location and channel, so the lambda function takes two arguments.
156 tbData = np.fromfunction(lambda iloc, ich: 250 + (((iloc - midLoc)**2 + (ich-midChan)**2))**0.5, (numLocs, numChans), dtype=float)
157 
158 # Write the data into the variables.
159 lonVar.writeNPArray.float(lonData)
160 latVar.writeNPArray.float(latData)
161 tbVar.writeNPArray.float(tbData)
static ObsGroup generate(Group &emptyGroup, const NewDimensionScales_t &fundamentalDims, std::shared_ptr< const detail::DataLayoutPolicy > layout=nullptr)
Create an empty ObsGroup and populate it with the fundamental dimensions.
Definition: ObsGroup.cpp:72
IODA_DL Group createFile(const std::string &filename, BackendCreateModes mode, HDF5_Version_Range compat=defaultVersionRange())
Create a ioda::Group backed by an HDF5 file.
Definition: HH.cpp:120
Used to specify Variable creation-time properties.
Definition: Has_Variables.h:57