IODA
obsspace_f.cc
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2017-2019 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 #include "ioda/core/obsspace_f.h"
9 
10 #include <algorithm>
11 #include <cstring>
12 #include <string>
13 #include <vector>
14 
15 #include "eckit/exception/Exceptions.h"
16 
17 #include "oops/util/DateTime.h"
18 
19 #include "ioda/ObsSpace.h"
20 
21 namespace ioda {
22 
23 // -----------------------------------------------------------------------------
24 const ObsSpace * obsspace_construct_f(const eckit::Configuration * conf,
25  const util::DateTime * begin,
26  const util::DateTime * end) {
27  return new ObsSpace(*conf, oops::mpi::world(), *begin, *end, oops::mpi::myself());
28 }
29 
30 // -----------------------------------------------------------------------------
32  ASSERT(obss != nullptr);
33  delete obss;
34 }
35 
36 // -----------------------------------------------------------------------------
37 void obsspace_obsname_f(const ObsSpace & obss, size_t & lcname, char * cname) {
38  std::string obsname = obss.obsname();
39  lcname = obsname.size();
40  ASSERT(lcname < 100); // to not overflow the associated fortran string
41  strncpy(cname, obsname.c_str(), lcname);
42 }
43 
44 // -----------------------------------------------------------------------------
45 const oops::Variables * obsspace_obsvariables_f(const ObsSpace & obss) {
46  return &obss.obsvariables();
47 }
48 
49 // -----------------------------------------------------------------------------
50 std::size_t obsspace_get_gnlocs_f(const ObsSpace & obss) {
51  return obss.globalNumLocs();
52 }
53 // -----------------------------------------------------------------------------
54 std::size_t obsspace_get_nlocs_f(const ObsSpace & obss) {
55  return obss.nlocs();
56 }
57 // -----------------------------------------------------------------------------
58 std::size_t obsspace_get_nchans_f(const ObsSpace & obss) {
59  return obss.nchans();
60 }
61 // -----------------------------------------------------------------------------
62 std::size_t obsspace_get_nrecs_f(const ObsSpace & obss) {
63  return obss.nrecs();
64 }
65 // -----------------------------------------------------------------------------
66 std::size_t obsspace_get_nvars_f(const ObsSpace & obss) {
67  return obss.nvars();
68 }
69 
70 // -----------------------------------------------------------------------------
71 void obsspace_get_dim_name_f(const ObsSpace & obss, const int & dim_id,
72  std::size_t & len_dim_name, char * dim_name) {
73  std::string dimName = obss.get_dim_name(static_cast<ioda::ObsDimensionId>(dim_id));
74  len_dim_name = dimName.size();
75  ASSERT(len_dim_name < 100); // to not overflow the associated fortran string
76  strncpy(dim_name, dimName.c_str(), len_dim_name);
77 }
78 // -----------------------------------------------------------------------------
79 std::size_t obsspace_get_dim_size_f(const ObsSpace & obss, const int & dim_id) {
80  return obss.get_dim_size(static_cast<ioda::ObsDimensionId>(dim_id));
81 }
82 // -----------------------------------------------------------------------------
83 int obsspace_get_dim_id_f(const ObsSpace & obss, const char * dim_name) {
84  return static_cast<int>(obss.get_dim_id(std::string(dim_name)));
85 }
86 
87 // -----------------------------------------------------------------------------
88 void obsspace_get_comm_f(const ObsSpace & obss, int & lcname, char * cname) {
89  lcname = obss.comm().name().size();
90  ASSERT(lcname < 100); // to not overflow the associated fortran string
91  strncpy(cname, obss.comm().name().c_str(), lcname);
92 }
93 // -----------------------------------------------------------------------------
94 void obsspace_get_recnum_f(const ObsSpace & obss,
95  const std::size_t & length, std::size_t * recnum) {
96  ASSERT(length >= obss.nlocs());
97  for (std::size_t i = 0; i < length; i++) {
98  recnum[i] = obss.recnum()[i];
99  }
100 }
101 // -----------------------------------------------------------------------------
102 void obsspace_get_index_f(const ObsSpace & obss,
103  const std::size_t & length, std::size_t * index) {
104  ASSERT(length >= obss.nlocs());
105  for (std::size_t i = 0; i < length; i++) {
106  // Fortran array indices start at 1, whereas C indices start at 0.
107  // Add 1 to each index value as it is handed off from C to Fortran.
108  index[i] = obss.index()[i] + 1;
109  }
110 }
111 // -----------------------------------------------------------------------------
112 bool obsspace_has_f(const ObsSpace & obss, const char * group, const char * vname) {
113  return obss.has(std::string(group), std::string(vname));
114 }
115 // -----------------------------------------------------------------------------
116 void obsspace_get_int32_f(const ObsSpace & obss, const char * group, const char * vname,
117  const std::size_t & length, int32_t* vec,
118  const std::size_t & len_cs, int* chan_select) {
119  ASSERT(len_cs <= obss.nchans());
120  std::vector<int> chanSelect(len_cs);
121  chanSelect.assign(chan_select, chan_select + len_cs);
122  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
123  else
124  ASSERT(length >= obss.nlocs());
125  std::vector<int32_t> vdata(length);
126  obss.get_db(std::string(group), std::string(vname), vdata, chanSelect);
127  std::copy(vdata.begin(), vdata.end(), vec);
128 }
129 // -----------------------------------------------------------------------------
130 void obsspace_get_int64_f(const ObsSpace & obss, const char * group, const char * vname,
131  const std::size_t & length, int64_t* vec,
132  const std::size_t & len_cs, int* chan_select) {
133  ASSERT(len_cs <= obss.nchans());
134  std::vector<int> chanSelect(len_cs);
135  chanSelect.assign(chan_select, chan_select + len_cs);
136  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
137  else
138  ASSERT(length >= obss.nlocs());
139 // obss.get_db(std::string(group), std::string(vname), vec, chanSelect);
140 }
141 // -----------------------------------------------------------------------------
142 void obsspace_get_real32_f(const ObsSpace & obss, const char * group, const char * vname,
143  const std::size_t & length, float* vec,
144  const std::size_t & len_cs, int* chan_select) {
145  ASSERT(len_cs <= obss.nchans());
146  std::vector<int> chanSelect(len_cs);
147  chanSelect.assign(chan_select, chan_select + len_cs);
148  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
149  else
150  ASSERT(length >= obss.nlocs());
151 // obss.get_db(std::string(group), std::string(vname), vec, chanSelect);
152 }
153 // -----------------------------------------------------------------------------
154 void obsspace_get_real64_f(const ObsSpace & obss, const char * group, const char * vname,
155  const std::size_t & length, double* vec,
156  const std::size_t & len_cs, int* chan_select) {
157  ASSERT(len_cs <= obss.nchans());
158  std::vector<int> chanSelect(len_cs);
159  chanSelect.assign(chan_select, chan_select + len_cs);
160  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
161  else
162  ASSERT(length >= obss.nlocs());
163  std::vector<double> vdata(length);
164  obss.get_db(std::string(group), std::string(vname), vdata, chanSelect);
165  std::copy(vdata.begin(), vdata.end(), vec);
166 }
167 // -----------------------------------------------------------------------------
168 void obsspace_get_datetime_f(const ObsSpace & obss, const char * group, const char * vname,
169  const std::size_t & length, int32_t* date, int32_t* time,
170  const std::size_t & len_cs, int* chan_select) {
171  ASSERT(len_cs <= obss.nchans());
172  std::vector<int> chanSelect(len_cs);
173  chanSelect.assign(chan_select, chan_select + len_cs);
174  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
175  else
176  ASSERT(length >= obss.nlocs());
177 
178  // Load a DateTime vector from the database, then convert to a date and time
179  // vector which are then returned.
180  util::DateTime temp_dt("0000-01-01T00:00:00Z");
181  std::vector<util::DateTime> dt_vect(length, temp_dt);
182  obss.get_db(std::string(group), std::string(vname), dt_vect, chanSelect);
183 
184  // Convert to date and time values. The DateTime utilities can return year, month,
185  // day, hour, minute second.
186  int year;
187  int month;
188  int day;
189  int hour;
190  int minute;
191  int second;
192  for (std::size_t i = 0; i < length; i++) {
193  dt_vect[i].toYYYYMMDDhhmmss(year, month, day, hour, minute, second);
194  date[i] = (year * 10000) + (month * 100) + day;
195  time[i] = (hour * 10000) + (minute * 100) + second;
196  }
197 }
198 // -----------------------------------------------------------------------------
199 void obsspace_put_int32_f(ObsSpace & obss, const char * group, const char * vname,
200  const std::size_t & length, int32_t* vec,
201  const std::size_t & ndims, int* dim_ids) {
202  // Generate the dimension list (names) and the total number of elements required
203  // (product of dimension sizes). vec is just an allocated memory buffer in which
204  // to place the data from the ObsSpace variable.
205  std::vector<std::string> dimList;
206  int numElements = (ndims > 0) ? 1 : 0;
207  for (std::size_t i = 0; i < ndims; ++i) {
208  ObsDimensionId dimId = static_cast<ioda::ObsDimensionId>(dim_ids[i]);
209  dimList.push_back(obss.get_dim_name(dimId));
210  numElements *= obss.get_dim_size(dimId);
211  }
212 
213  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
214  else
215  ASSERT(length >= numElements);
216  std::vector<int32_t> vdata;
217  vdata.assign(vec, vec + length);
218 
219  obss.put_db(std::string(group), std::string(vname), vdata, dimList);
220 }
221 // -----------------------------------------------------------------------------
222 void obsspace_put_int64_f(ObsSpace & obss, const char * group, const char * vname,
223  const std::size_t & length, int64_t* vec,
224  const std::size_t & ndims, int* dim_ids) {
225  // Generate the dimension list (names) and the total number of elements required
226  // (product of dimension sizes). vec is just an allocated memory buffer in which
227  // to place the data from the ObsSpace variable.
228  std::vector<std::string> dimList;
229  int numElements = (ndims > 0) ? 1 : 0;
230  for (std::size_t i = 0; i < ndims; ++i) {
231  ObsDimensionId dimId = static_cast<ioda::ObsDimensionId>(dim_ids[i]);
232  dimList.push_back(obss.get_dim_name(dimId));
233  numElements *= obss.get_dim_size(dimId);
234  }
235 
236  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
237  else
238  ASSERT(length >= numElements);
239 // obss.put_db(std::string(group), std::string(vname), vec, dimList);
240 }
241 // -----------------------------------------------------------------------------
242 void obsspace_put_real32_f(ObsSpace & obss, const char * group, const char * vname,
243  const std::size_t & length, float* vec,
244  const std::size_t & ndims, int* dim_ids) {
245  // Generate the dimension list (names) and the total number of elements required
246  // (product of dimension sizes). vec is just an allocated memory buffer in which
247  // to place the data from the ObsSpace variable.
248  std::vector<std::string> dimList;
249  int numElements = (ndims > 0) ? 1 : 0;
250  for (std::size_t i = 0; i < ndims; ++i) {
251  ObsDimensionId dimId = static_cast<ioda::ObsDimensionId>(dim_ids[i]);
252  dimList.push_back(obss.get_dim_name(dimId));
253  numElements *= obss.get_dim_size(dimId);
254  }
255 
256  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
257  else
258  ASSERT(length >= numElements);
259 // obss.put_db(std::string(group), std::string(vname), vec, dimList);
260 }
261 // -----------------------------------------------------------------------------
262 void obsspace_put_real64_f(ObsSpace & obss, const char * group, const char * vname,
263  const std::size_t & length, double* vec,
264  const std::size_t & ndims, int* dim_ids) {
265  // Generate the dimension list (names) and the total number of elements required
266  // (product of dimension sizes). vec is just an allocated memory buffer in which
267  // to place the data from the ObsSpace variable.
268  std::vector<std::string> dimList;
269  int numElements = (ndims > 0) ? 1 : 0;
270  for (std::size_t i = 0; i < ndims; ++i) {
271  ObsDimensionId dimId = static_cast<ioda::ObsDimensionId>(dim_ids[i]);
272  dimList.push_back(obss.get_dim_name(dimId));
273  numElements *= obss.get_dim_size(dimId);
274  }
275 
276  if (std::string(group) == "VarMetaData") ASSERT(length >= obss.nvars());
277  else
278  ASSERT(length >= numElements);
279  std::vector<double> vdata;
280  vdata.assign(vec, vec + length);
281 
282  obss.put_db(std::string(group), std::string(vname), vdata, dimList);
283 }
284 // -----------------------------------------------------------------------------
286  return static_cast<int>(ObsDimensionId::Nlocs);
287 }
288 // -----------------------------------------------------------------------------
290  return static_cast<int>(ObsDimensionId::Nchans);
291 }
292 // -----------------------------------------------------------------------------
293 
294 } // namespace ioda
Observation data class for IODA.
Definition: src/ObsSpace.h:116
std::size_t get_dim_size(const ObsDimensionId dimId) const
return the standard dimension size for the given dimension id
Definition: src/ObsSpace.h:197
bool has(const std::string &group, const std::string &name) const
return true if group/variable exists
Definition: ObsSpace.cc:230
const std::vector< std::size_t > & recnum() const
return reference to the record number vector
Definition: src/ObsSpace.h:222
void put_db(const std::string &group, const std::string &name, const std::vector< int > &vdata, const std::vector< std::string > &dimList={ "nlocs" })
transfer data from vdata to the obs container
Definition: ObsSpace.cc:306
size_t nchans() const
return the number of channels in the container. If this is not a radiance obs type,...
Definition: src/ObsSpace.h:180
std::string get_dim_name(const ObsDimensionId dimId) const
return the standard dimension name for the given dimension id
Definition: src/ObsSpace.h:192
std::size_t nrecs() const
return the number of records in the obs space container
Definition: src/ObsSpace.h:185
const std::vector< std::size_t > & index() const
return reference to the index vector
Definition: src/ObsSpace.h:241
const eckit::mpi::Comm & comm() const
Definition: src/ObsSpace.h:148
ObsDimensionId get_dim_id(const std::string &dimName) const
return the standard dimension id for the given dimension name
Definition: src/ObsSpace.h:202
size_t nlocs() const
return the number of locations in the obs space. Note that nlocs may be smaller than global unique nl...
Definition: src/ObsSpace.h:176
void get_db(const std::string &group, const std::string &name, std::vector< int > &vdata, const std::vector< int > &chanSelect={ }) const
transfer data from the obs container to vdata
Definition: ObsSpace.cc:270
std::size_t nvars() const
return the number of variables in the obs space container. "Variables" refers to the quantities that ...
Definition: ObsSpace.cc:189
std::size_t globalNumLocs() const
return the total number of locations in the corresponding obs spaces across all MPI tasks
Definition: src/ObsSpace.h:168
const oops::Variables & obsvariables() const
return oops variables object (simulated variables)
Definition: src/ObsSpace.h:332
const std::string & obsname() const
return the name of the obs type being stored
Definition: src/ObsSpace.h:216
std::size_t obsspace_get_nvars_f(const ObsSpace &obss)
Definition: obsspace_f.cc:66
void obsspace_get_comm_f(const ObsSpace &obss, int &lcname, char *cname)
Definition: obsspace_f.cc:88
const ObsSpace * obsspace_construct_f(const eckit::Configuration *conf, const util::DateTime *begin, const util::DateTime *end)
Definition: obsspace_f.cc:24
void obsspace_get_int32_f(const ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, int32_t *vec, const std::size_t &len_cs, int *chan_select)
Definition: obsspace_f.cc:116
void obsspace_put_real32_f(ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, float *vec, const std::size_t &ndims, int *dim_ids)
Definition: obsspace_f.cc:242
void obsspace_put_int32_f(ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, int32_t *vec, const std::size_t &ndims, int *dim_ids)
Definition: obsspace_f.cc:199
void obsspace_get_real32_f(const ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, float *vec, const std::size_t &len_cs, int *chan_select)
Definition: obsspace_f.cc:142
void obsspace_destruct_f(ObsSpace *obss)
Definition: obsspace_f.cc:31
std::size_t obsspace_get_nlocs_f(const ObsSpace &obss)
Definition: obsspace_f.cc:54
const oops::Variables * obsspace_obsvariables_f(const ObsSpace &obss)
Definition: obsspace_f.cc:45
IODA_DL void copy(const ObjectSelection &from, ObjectSelection &to, const ScaleMapping &scale_map)
Generic data copying function.
Definition: Copying.cpp:63
void obsspace_get_recnum_f(const ObsSpace &obss, const std::size_t &length, std::size_t *recnum)
Definition: obsspace_f.cc:94
int obsspace_get_dim_id_f(const ObsSpace &obss, const char *dim_name)
Definition: obsspace_f.cc:83
void obsspace_obsname_f(const ObsSpace &obss, size_t &lcname, char *cname)
Definition: obsspace_f.cc:37
void obsspace_get_real64_f(const ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, double *vec, const std::size_t &len_cs, int *chan_select)
Definition: obsspace_f.cc:154
void obsspace_get_datetime_f(const ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, int32_t *date, int32_t *time, const std::size_t &len_cs, int *chan_select)
Definition: obsspace_f.cc:168
void obsspace_get_index_f(const ObsSpace &obss, const std::size_t &length, std::size_t *index)
Definition: obsspace_f.cc:102
void obsspace_put_real64_f(ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, double *vec, const std::size_t &ndims, int *dim_ids)
Definition: obsspace_f.cc:262
bool obsspace_has_f(const ObsSpace &obss, const char *group, const char *vname)
Definition: obsspace_f.cc:112
void obsspace_get_int64_f(const ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, int64_t *vec, const std::size_t &len_cs, int *chan_select)
Definition: obsspace_f.cc:130
std::size_t obsspace_get_nchans_f(const ObsSpace &obss)
Definition: obsspace_f.cc:58
int obsspace_get_nchans_dim_id_f()
Definition: obsspace_f.cc:289
int obsspace_get_nlocs_dim_id_f()
Definition: obsspace_f.cc:285
void obsspace_put_int64_f(ObsSpace &obss, const char *group, const char *vname, const std::size_t &length, int64_t *vec, const std::size_t &ndims, int *dim_ids)
Definition: obsspace_f.cc:222
std::size_t obsspace_get_dim_size_f(const ObsSpace &obss, const int &dim_id)
Definition: obsspace_f.cc:79
std::size_t obsspace_get_gnlocs_f(const ObsSpace &obss)
Definition: obsspace_f.cc:50
void obsspace_get_dim_name_f(const ObsSpace &obss, const int &dim_id, std::size_t &len_dim_name, char *dim_name)
Definition: obsspace_f.cc:71
ObsDimensionId
Definition: src/ObsSpace.h:61
std::size_t obsspace_get_nrecs_f(const ObsSpace &obss)
Definition: obsspace_f.cc:62