TECA
The Toolkit for Extreme Climate Analysis
teca_coordinate_util.h
Go to the documentation of this file.
1 #ifndef teca_cartesian_mesh_util_h
2 #define teca_cartesian_mesh_util_h
3 
4 /// @file
5 
6 #include "teca_config.h"
7 #include "teca_cartesian_mesh.h"
8 #include "teca_variant_array.h"
10 #include "teca_metadata.h"
11 #include "teca_array_attributes.h"
12 
13 #include <vector>
14 #include <cmath>
15 #include <type_traits>
16 #include <typeinfo>
17 #include <iomanip>
18 
19 #if defined(TECA_HAS_CUDA)
20 #define TECA_TARGET __host__ __device__
21 #else
22 #define TECA_TARGET
23 #endif
24 
25 /// For printing data as ASCII with the maximum supported numerical precision
26 #define max_prec(T) \
27  std::setprecision(std::numeric_limits<T>::digits10 + 1)
28 
29 /// Codes dealing with operations on coordinate systems
31 {
32 /** @brief
33  * traits classes used to get default tolerances for comparing numbers
34  * of a given precision.
35  *
36  * @details
37  * A relative tolerance is used for comparing large
38  * numbers and an absolute tolerance is used for comparing small numbers.
39  * these defaults are not universal and will not work well in all situations.
40  *
41  * see also:
42  * https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
43  */
44 template <typename n_t>
45 struct equal_tt {};
46 
47 #define declare_equal_tt(cpp_t, atol, rtol) \
48 /** Specialization for cpp_t with default absTol and relTol */ \
49 template <> \
50 struct equal_tt<cpp_t> \
51 { \
52  TECA_TARGET static cpp_t absTol() { return atol; } \
53  TECA_TARGET static cpp_t relTol() { return rtol; } \
54 };
55 
56 declare_equal_tt(float, 10.0f*std::numeric_limits<float>::epsilon(),
57  std::numeric_limits<float>::epsilon())
58 
59 declare_equal_tt(double, 10.0*std::numeric_limits<double>::epsilon(),
60  std::numeric_limits<float>::epsilon())
61 
62 declare_equal_tt(long double, std::numeric_limits<double>::epsilon(),
63  std::numeric_limits<double>::epsilon())
64 
65 /** Compare two floating point numbers. absTol handles comparing numbers very
66  * close to zero. relTol handles comparing larger values.
67  */
68 template <typename T>
69 TECA_TARGET
70 bool equal(T a, T b,
71  T relTol = equal_tt<T>::relTol(), T absTol = equal_tt<T>::absTol(),
72  typename std::enable_if<std::is_floating_point<T>::value>::type* = 0)
73 {
74  // for numbers close to zero
75  T diff = std::abs(a - b);
76  if (diff <= absTol)
77  return true;
78  // realtive difference for larger values
79  a = std::abs(a);
80  b = std::abs(b);
81  b = (b > a) ? b : a;
82  b *= relTol;
83  if (diff <= b)
84  return true;
85  return false;
86 }
87 
88 /// Compare two integral numbers.
89 template <typename T>
90 TECA_TARGET
91 bool equal(T a, T b, T relTol = 0, T absTol = 0,
92  typename std::enable_if<std::is_integral<T>::value>::type* = 0)
93 {
94  (void)relTol;
95  (void)absTol;
96  return a == b;
97 }
98 
99 /** Compare two floating point numbers. This overload may be used in regression
100  * tests or other contexts where a diagnostic error message should be reported
101  * if the numbers are not equal.
102  */
103 template <typename T>
104 bool equal(T a, T b, std::string &diagnostic,
105  T relTol = equal_tt<T>::relTol(), T absTol = equal_tt<T>::absTol(),
106  typename std::enable_if<std::is_floating_point<T>::value>::type* = 0)
107 {
108  // for numbers close to zero
109  T diff = std::abs(a - b);
110  if (diff <= absTol)
111  return true;
112  // realtive difference for larger values
113  T aa = std::abs(a);
114  T bb = std::abs(b);
115  bb = (bb > aa) ? bb : aa;
116  T cc = bb*relTol;
117  if (diff <= cc)
118  return true;
119  // a and b are not equal format the diagnostic
120  T eps = std::numeric_limits<T>::epsilon();
121  std::ostringstream os;
122  os << max_prec(T) << a << " != " << max_prec(T) << b
123  << " relTol=" << max_prec(T) << relTol
124  << " absTol=" << max_prec(T) << absTol
125  << " |a-b|=" << max_prec(T) << diff
126  << " |a-b|/eps=" << max_prec(T) << diff/eps
127  << " max(|a|,|b|)*relTol=" << max_prec(T) << cc
128  << " |a-b|/max(|a|,|b|)=" << max_prec(T) << diff/bb
129  << " eps(" << typeid(a).name() << sizeof(T) << ")="
130  << max_prec(T) << eps;
131  diagnostic = os.str();
132  return false;
133 }
134 
135 /** Compare two integral numbers. This overload may be used in regression
136  * tests or other contexts where a diagnostic error message should be reported
137  * if the numbers are not equal.
138  */
139 template <typename T>
140 bool equal(T a, T b, std::string &diagnostic, T relTol = 0, T absTol = 0,
141  typename std::enable_if<std::is_integral<T>::value>::type* = 0)
142 {
143  (void)relTol;
144  (void)absTol;
145  if (a == b)
146  return true;
147  // a and b are not equal format the diagnostic
148  std::ostringstream os;
149  os << typeid(a).name() << sizeof(T) << " " << a << " != " << b;
150  diagnostic = os.str();
151  return false;
152 }
153 
154 /** Error codes returned by bool equal(const const_p_teca_variant_array &array1,
155  const const_p_teca_variant_array &array2,
156  double absTol, double relTol, int errorNo,
157  std::string &errorStr);
158  * | Code | Meaning |
159  * |------------------|---------|
160  * | no_error | no error |
161  * | length_missmatch | length missmatch |
162  * | type_missmatch | type missmatch |
163  * | invalid_value | invalid value detected (NaN, inf) |
164  * | value_missmatch | a value failed to compare within the tolerance |
165  */
167 {
168  enum {no_error = 0,
169  invalid_value = 1,
170  length_missmatch = 2,
171  type_missmatch = 3,
172  value_missmatch = 4,
173  unsupported_type = 5
174  };
175 };
176 
177 /** Compare two variant arrays elementwise for equality. If the arrays fail to
178  * compare within the specified tolerance errorNo will contain one of the
179  * equal_error enumerations and errorStr will conatin a diagnostic message
180  * describing the failure.
181  */
183 bool equal(const const_p_teca_variant_array &array1,
184  const const_p_teca_variant_array &array2,
185  double absTol, double relTol, int &errorNo,
186  std::string &errorStr);
187 
188 /// Less than or equal to predicate
189 template<typename data_t>
191 { static bool eval(const data_t &l, const data_t &r) { return l <= r; } };
192 
193 /// Greater than or equal to predicate
194 template<typename data_t>
196 { static bool eval(const data_t &l, const data_t &r) { return l >= r; } };
197 
198 /// Less than predicate
199 template<typename data_t>
201 { static bool eval(const data_t &l, const data_t &r) { return l < r; } };
202 
203 /// Greater than predicate
204 template<typename data_t>
206 { static bool eval(const data_t &l, const data_t &r) { return l > r; } };
207 
208 /// comparator implementing bracket for ascending input arrays
209 template<typename data_t>
211 {
212  // m_0 is an index into the data, m_1 = m_0 + 1
213  // comparitors defining the bracket orientation. for data in
214  // ascending order: val >= data[m_0] && val <= data[m_1]
215  using comp0_t = geq<data_t>;
216  using comp1_t = leq<data_t>;
217 
218  // m_0 is an index into the data, m_1 = m_0 + 1
219  // get the id of the smaller value (lower == true)
220  // or the larger value (lower == false)
221  static unsigned long get_id(bool lower,
222  unsigned long m_0, unsigned long m_1)
223  {
224  if (lower)
225  return m_0;
226  return m_1;
227  }
228 };
229 
230 /// comparator implementing bracket for descending input arrays
231 template<typename data_t>
233 {
234  // m_0 is an index into the data, m_1 = m_0 + 1
235  // comparitors defining the bracket orientation. for data in
236  // descending order: val <= data[m_0] && val >= data[m_1]
237  using comp0_t = leq<data_t>;
238  using comp1_t = geq<data_t>;
239 
240  // m_0 is an index into the data, m_1 = m_0 + 1
241  // get the id of the smaller value (lower == true)
242  // or the larger value (lower == false)
243  static unsigned long get_id(bool lower,
244  unsigned long m_0, unsigned long m_1)
245  {
246  if (lower)
247  return m_1;
248  return m_0;
249  }
250 };
251 
252 /** binary search that will locate index bounding the value above or below
253  * such that data[i] <= val or val <= data[i+1] depending on the value of
254  * lower. return 0 if the value is found. the comp0 and comp1 template
255  * parameters let us operate on both ascending and descending input. defaults
256  * are set for ascending inputs.
257  */
258 template <typename data_t, typename bracket_t = ascend_bracket<data_t>>
260 int index_of(const data_t *data, unsigned long l, unsigned long r,
261  data_t val, bool lower, unsigned long &id)
262 {
263  unsigned long m_0 = (r + l)/2;
264  unsigned long m_1 = m_0 + 1;
265 
266  if (m_0 == r)
267  {
268  if (equal(val, data[m_0]))
269  {
270  id = m_0;
271  return 0;
272  }
273  // not found
274  return -1;
275  }
276  else
277  if (bracket_t::comp0_t::eval(val, data[m_0]) &&
278  bracket_t::comp1_t::eval(val, data[m_1]))
279  {
280  // found a bracket around the value
281  if (equal(val, data[m_0]))
282  id = m_0;
283  else
284  if (equal(val, data[m_1]))
285  id = m_1;
286  else
287  id = bracket_t::get_id(lower, m_0, m_1);
288  return 0;
289  }
290  else
291  if (bracket_t::comp1_t::eval(val, data[m_0]))
292  {
293  // split range to the left
294  return teca_coordinate_util::index_of<data_t, bracket_t>(
295  data, l, m_0, val, lower, id);
296  }
297  else
298  {
299  // split the range to the right
300  return teca_coordinate_util::index_of<data_t, bracket_t>(
301  data, m_1, r, val, lower, id);
302  }
303 
304  // not found
305  return -1;
306 }
307 
308 /** binary search that will locate index of the given value. return 0 if the
309  * value is found.
310  */
311 template <typename T>
313 int index_of(const T *data, size_t l, size_t r, T val, unsigned long &id)
314 {
315  unsigned long m_0 = (r + l)/2;
316  unsigned long m_1 = m_0 + 1;
317 
318  if (m_0 == r)
319  {
320  // need this test when len of data is 1
321  if (equal(val, data[m_0], T(8)*std::numeric_limits<T>::epsilon()))
322  {
323  id = m_0;
324  return 0;
325  }
326  // not found
327  return -1;
328  }
329  else
330  if (equal(val, data[m_0], T(8)*std::numeric_limits<T>::epsilon()))
331  {
332  id = m_0;
333  return 0;
334  }
335  else
336  if (equal(val, data[m_1], T(8)*std::numeric_limits<T>::epsilon()))
337  {
338  id = m_1;
339  return 0;
340  }
341  else
342  if (val < data[m_0])
343  {
344  // split range to the left
345  return teca_coordinate_util::index_of(data, l, m_0, val, id);
346  }
347  else
348  {
349  // split the range to the right
350  return teca_coordinate_util::index_of(data, m_1, r, val, id);
351  }
352 
353  // not found
354  return -1;
355 }
356 
357 /** Convert bounds to extents. return non-zero if the requested bounds are
358  * not in the given coordinate arrays. coordinate arrays must not be empty.
359  */
361 int bounds_to_extent(const double *bounds,
363  const const_p_teca_variant_array &z, unsigned long *extent);
364 
366 int bounds_to_extent(const double *bounds,
367  const const_p_teca_variant_array &x, unsigned long *extent);
368 
370 int bounds_to_extent(const double *bounds, const teca_metadata &md,
371  unsigned long *extent);
372 
373 /** Get the i,j,k cell index of point x,y,z in the given mesh. return 0 if
374  * successful.
375  */
376 template<typename T>
378 int index_of(const const_p_teca_cartesian_mesh &mesh, T x, T y, T z,
379  unsigned long &i, unsigned long &j, unsigned long &k)
380 {
381  const_p_teca_variant_array xc = mesh->get_x_coordinates();
382  const_p_teca_variant_array yc = mesh->get_y_coordinates();
383  const_p_teca_variant_array zc = mesh->get_z_coordinates();
384 
387  xc.get(),
388 
389  auto p_xc = std::dynamic_pointer_cast<TT>(xc)->get_cpu_accessible();
390  auto p_yc = std::dynamic_pointer_cast<TT>(yc)->get_cpu_accessible();
391  auto p_zc = std::dynamic_pointer_cast<TT>(zc)->get_cpu_accessible();
392 
393  unsigned long nx = xc->size();
394  unsigned long ny = yc->size();
395  unsigned long nz = zc->size();
396 
397  if (teca_coordinate_util::index_of(p_xc.get(), 0, nx-1, x, true, i)
398  || teca_coordinate_util::index_of(p_yc.get(), 0, ny-1, y, true, j)
399  || teca_coordinate_util::index_of(p_zc.get(), 0, nz-1, z, true, k))
400  {
401  // out of bounds
402  return -1;
403  }
404 
405  // success
406  return 0;
407  )
408 
409  // coords are not a floating point type
410  return -1;
411 }
412 
413 /** given a human readable date string in YYYY-MM-DD hh:mm:ss format and a
414  * list of floating point offset times in the specified calendar and units find
415  * the closest time step. return 0 if successful see index_of for a description
416  * of lower, if clamp is true then when the date falls outside of the time
417  * values either the first or last time step is returned.
418  */
421  bool lower, bool clamp, const std::string &calendar,
422  const std::string &units, const std::string &date,
423  unsigned long &step);
424 
425 /** given a time value (val), associated time units (units), and calendar
426  * (calendar), return a human-readable rendering of the date (date) in a
427  * strftime-format (format). return 0 if successful.
428  */
430 int time_to_string(double val, const std::string &calendar,
431  const std::string &units, const std::string &format, std::string &date);
432 
433 /** build random access data structures for an indexed table. the index column
434  * gives each entity a unique id. the index is used to identify rows that
435  * belong in the entity. it is assumed that an entity occupies consecutive rows.
436  * the returns are: n_entities, the number of entities found; counts, the
437  * number of rows used by each entity; offsets, the starting row of each
438  * entity; ids, a new set of ids for the entities starting from 0
439  */
440 template <typename int_t>
442 void get_table_offsets(const int_t *index, unsigned long n_rows,
443  unsigned long &n_entities, std::vector<unsigned long> &counts,
444  std::vector<unsigned long> &offsets, std::vector<unsigned long> &ids)
445 {
446  // count unique number of steps and compute an array index
447  // from each time step.
448  n_entities = 1;
449  ids.resize(n_rows);
450  unsigned long n_m1 = n_rows - 1;
451  for (unsigned long i = 0; i < n_m1; ++i)
452  {
453  ids[i] = n_entities - 1;
454  if (index[i] != index[i+1])
455  ++n_entities;
456  }
457  ids[n_m1] = n_entities - 1;
458 
459  // compute num storms in each step
460  counts.resize(n_entities);
461  unsigned long q = 0;
462  for (unsigned long i = 0; i < n_entities; ++i)
463  {
464  counts[i] = 1;
465  while ((q < n_m1) && (index[q] == index[q+1]))
466  {
467  ++counts[i];
468  ++q;
469  }
470  ++q;
471  }
472 
473  // compute the offset to the first storm in each step
474  offsets.resize(n_entities);
475  offsets[0] = 0;
476  for (unsigned long i = 1; i < n_entities; ++i)
477  offsets[i] = offsets[i-1] + counts[i-1];
478 }
479 
480 /** 0th order (nearest neighbor) interpolation for nodal data on a stretched
481  * Cartesian mesh. This overload implements the general 3D case.
482  * cx, cy, cz is the location to interpolate to
483  * p_x, p_y, p_z array arrays containing the source coordinates with extents
484  * [0, ihi, 0, jhi, 0, khi]
485  * p_data is the field to interpolate from
486  * val is the result
487  * returns 0 if successful, an error occurs if cx, cy, cz is outside of the
488  * source coordinate system
489  */
490 template<typename CT, typename DT>
492 int interpolate_nearest(CT cx, CT cy, CT cz,
493  const CT *p_x, const CT *p_y, const CT *p_z,
494  const DT *p_data, unsigned long ihi, unsigned long jhi,
495  unsigned long khi, unsigned long nx, unsigned long nxy,
496  DT &val)
497 {
498  // get i,j of node less than cx,cy
499  unsigned long i = 0;
500  unsigned long j = 0;
501  unsigned long k = 0;
502 
503  if ((ihi && teca_coordinate_util::index_of(p_x, 0, ihi, cx, true, i))
504  || (jhi && teca_coordinate_util::index_of(p_y, 0, jhi, cy, true, j))
505  || (khi && teca_coordinate_util::index_of(p_z, 0, khi, cz, true, k)))
506  {
507  // cx,cy,cz is outside the coordinate axes
508  return -1;
509  }
510 
511  // get i,j of node greater than cx,cy
512  unsigned long ii = std::min(i + 1, ihi);
513  unsigned long jj = std::min(j + 1, jhi);
514  unsigned long kk = std::min(k + 1, khi);
515 
516  // get index of nearest node
517  unsigned long p = (cx - p_x[i]) <= (p_x[ii] - cx) ? i : ii;
518  unsigned long q = (cy - p_y[j]) <= (p_y[jj] - cy) ? j : jj;
519  unsigned long r = (cz - p_z[k]) <= (p_z[kk] - cz) ? k : kk;
520 
521  // assign value from nearest node
522  val = p_data[p + nx*q + nxy*r];
523  return 0;
524 }
525 
526 /** 0th order (nearest neighbor) interpolation for nodal data on a stretched
527  * Cartesian mesh. This overload implements the special case where both source
528  * and target mesh data are in a 2D x-y plane using fewer operations than the
529  * general 3D implementation.
530  * cx, cy, cz is the location to interpolate to
531  * p_x, p_y, p_z array arrays containing the source coordinates with extents
532  * [0, ihi, 0, jhi, 0, khi]
533  * p_data is the field to interpolate from
534  * val is the result
535  * returns 0 if successful, an error occurs if cx, cy, cz is outside of the
536  * source coordinate system
537  */
538 template<typename coord_t, typename data_t>
540 int interpolate_nearest(coord_t cx, coord_t cy, const coord_t *p_x,
541  const coord_t *p_y, const data_t *p_data, unsigned long ihi,
542  unsigned long jhi, unsigned long nx, data_t &val)
543 {
544  // get i,j of node less than cx,cy
545  unsigned long i = 0;
546  unsigned long j = 0;
547 
548  if ((ihi && teca_coordinate_util::index_of(p_x, 0, ihi, cx, true, i))
549  || (jhi && teca_coordinate_util::index_of(p_y, 0, jhi, cy, true, j)))
550  {
551  // cx,cy is outside the coordinate axes
552  return -1;
553  }
554 
555  // get i,j of node greater than cx,cy
556  unsigned long ii = std::min(i + 1, ihi);
557  unsigned long jj = std::min(j + 1, jhi);
558 
559  // get index of nearest node
560  unsigned long p = (cx - p_x[i]) <= (p_x[ii] - cx) ? i : ii;
561  unsigned long q = (cy - p_y[j]) <= (p_y[jj] - cy) ? j : jj;
562 
563  // assign value from nearest node
564  val = p_data[p + nx*q];
565 
566  return 0;
567 }
568 
569 /** 1st order (linear) interpolation for nodal data on stretched Cartesian
570  * mesh. This overload implements the general 3D case.
571  * cx, cy, cz is the location to interpolate to
572  * p_x, p_y, p_z array arrays containing the source coordinates with extents
573  * [0, ihi, 0, jhi, 0, khi]
574  * p_data is the field to interpolate from
575  * val is the result
576  * returns 0 if successful, an error occurs if cx, cy, cz is outside of the
577  * source coordinate system
578  */
579 template<typename CT, typename DT>
581 int interpolate_linear(CT cx, CT cy, CT cz,
582  const CT *p_x, const CT *p_y, const CT *p_z,
583  const DT *p_data, unsigned long ihi, unsigned long jhi,
584  unsigned long khi, unsigned long nx, unsigned long nxy,
585  DT &val)
586 {
587  // get i,j of node less than cx,cy
588  unsigned long i = 0;
589  unsigned long j = 0;
590  unsigned long k = 0;
591 
592  if ((ihi && teca_coordinate_util::index_of(p_x, 0, ihi, cx, true, i))
593  || (jhi && teca_coordinate_util::index_of(p_y, 0, jhi, cy, true, j))
594  || (khi && teca_coordinate_util::index_of(p_z, 0, khi, cz, true, k)))
595  {
596  // cx,cy,cz is outside the coordinate axes
597  return -1;
598  }
599 
600  // get i,j of node greater than cx,cy,cz
601  unsigned long ii = std::min(i + 1, ihi);
602  unsigned long jj = std::min(j + 1, jhi);
603  unsigned long kk = std::min(k + 1, khi);
604 
605  // compute weights
606  CT wx = ii == i ? 0 : (cx - p_x[i])/(p_x[ii] - p_x[i]);
607  CT wy = jj == j ? 0 : (cy - p_y[j])/(p_y[jj] - p_y[j]);
608  CT wz = kk == k ? 0 : (cz - p_z[k])/(p_z[kk] - p_z[k]);
609 
610  CT vx = CT(1) - wx;
611  CT vy = CT(1) - wy;
612  CT vz = CT(1) - wz;
613 
614  // interpolate
615  val = vx*vy*vz*p_data[ i + j*nx + k*nxy]
616  + wx*vy*vz*p_data[ii + j*nx + k*nxy]
617  + wx*wy*vz*p_data[ii + jj*nx + k*nxy]
618  + vx*wy*vz*p_data[ i + jj*nx + k*nxy]
619  + vx*vy*wz*p_data[ i + j*nx + kk*nxy]
620  + wx*vy*wz*p_data[ii + j*nx + kk*nxy]
621  + wx*wy*wz*p_data[ii + jj*nx + kk*nxy]
622  + vx*wy*wz*p_data[ i + jj*nx + kk*nxy];
623 
624  return 0;
625 }
626 
627 /** 1st order (linear) interpolation for nodal data on stretched Cartesian mesh.
628  * This overload implements the special case where both source and target data
629  * are in a 2D x-y plane using fewer operations than the general 3D
630  * implementation.
631  * cx, cy, cz is the location to interpolate to
632  * p_x, p_y, p_z array arrays containing the source coordinates with extents
633  * [0, ihi, 0, jhi, 0, khi]
634  * p_data is the field to interpolate from
635  * val is the result
636  * returns 0 if successful, an error occurs if cx, cy, cz is outside of the
637  * source coordinate system
638  */
639 template<typename CT, typename DT>
641 int interpolate_linear(CT cx, CT cy, const CT *p_x, const CT *p_y,
642  const DT *p_data, unsigned long ihi, unsigned long jhi,
643  unsigned long nx, DT &val)
644 {
645  // get i,j of node less than cx,cy
646  unsigned long i = 0;
647  unsigned long j = 0;
648 
649  if ((ihi && teca_coordinate_util::index_of(p_x, 0, ihi, cx, true, i))
650  || (jhi && teca_coordinate_util::index_of(p_y, 0, jhi, cy, true, j)))
651  {
652  // cx,cy is outside the coordinate axes
653  return -1;
654  }
655 
656  // get i,j of node greater than cx,cy
657  unsigned long ii = std::min(i + 1, ihi);
658  unsigned long jj = std::min(j + 1, jhi);
659 
660  // compute weights
661  CT wx = ii == i ? 0 : (cx - p_x[i])/(p_x[ii] - p_x[i]);
662  CT wy = jj == j ? 0 : (cy - p_y[j])/(p_y[jj] - p_y[j]);
663 
664  CT vx = CT(1) - wx;
665  CT vy = CT(1) - wy;
666 
667  // interpolate
668  val = vx*vy*p_data[ i + j*nx]
669  + wx*vy*p_data[ii + j*nx]
670  + wx*wy*p_data[ii + jj*nx]
671  + vx*wy*p_data[ i + jj*nx];
672 
673  return 0;
674 }
675 
676 /// A functor templated on order of accuracy for above Cartesian mesh interpolants
677 template<int> struct TECA_EXPORT interpolate_t;
678 
679 /// Zero'th order interpolant specialization
680 template<> struct interpolate_t<0>
681 {
682  // 3D
683  template<typename CT, typename DT>
684  int operator()(CT tx, CT ty, CT tz, const CT *sx, const CT *sy,
685  const CT *sz, const DT *sa, unsigned long ihi, unsigned long jhi,
686  unsigned long khi, unsigned long nx, unsigned long nxy, DT &ta)
687  {
689  sx,sy,sz,sa, ihi,jhi,khi, nx,nxy, ta);
690  }
691 
692  // 2D x-y plane
693  template<typename CT, typename DT>
694  int operator()(CT tx, CT ty, const CT *sx, const CT *sy,
695  const DT *sa, unsigned long ihi, unsigned long jhi,
696  unsigned long nx, DT &ta)
697  {
699  sx,sy,sa, ihi,jhi, nx, ta);
700  }
701 };
702 
703 /// First order interpolant specialization
704 template<> struct interpolate_t<1>
705 {
706  // 3D
707  template<typename CT, typename DT>
708  int operator()(CT tx, CT ty, CT tz, const CT *sx, const CT *sy,
709  const CT *sz, const DT *sa, unsigned long ihi, unsigned long jhi,
710  unsigned long khi, unsigned long nx,unsigned long nxy, DT &ta)
711  {
713  sx,sy,sz,sa, ihi,jhi,khi, nx,nxy, ta);
714  }
715 
716  // 2D x-y plane
717  template<typename CT, typename DT>
718  int operator()(CT tx, CT ty, const CT *sx, const CT *sy,
719  const DT *sa, unsigned long ihi, unsigned long jhi,
720  unsigned long nx, DT &ta)
721  {
723  sx,sy,sa, ihi,jhi, nx, ta);
724  }
725 };
726 
727 /// return 0 if the centering is one of the values defined in teca_array_attributes
729 int validate_centering(int centering);
730 
731 /// convert from a cell extent to a face, edge or point centered extent
732 template <typename num_t>
734 int convert_cell_extent(num_t *extent, int centering)
735 {
736  switch (centering)
737  {
738  case teca_array_attributes::invalid_value:
739  TECA_ERROR("detected invalid_value in centering")
740  return -1;
741  break;
742  case teca_array_attributes::cell_centering:
743  break;
744  case teca_array_attributes::x_face_centering:
745  extent[1] += 1;
746  break;
747  case teca_array_attributes::y_face_centering:
748  extent[3] += 1;
749  break;
750  case teca_array_attributes::z_face_centering:
751  extent[5] += 1;
752  break;
753  case teca_array_attributes::x_edge_centering:
754  extent[3] += 1;
755  extent[5] += 1;
756  break;
757  case teca_array_attributes::y_edge_centering:
758  extent[1] += 1;
759  extent[5] += 1;
760  break;
761  case teca_array_attributes::z_edge_centering:
762  extent[1] += 1;
763  extent[3] += 1;
764  break;
765  case teca_array_attributes::point_centering:
766  extent[1] += 1;
767  extent[3] += 1;
768  extent[5] += 1;
769  break;
770  case teca_array_attributes::no_centering:
771  break;
772  default:
773  TECA_ERROR("this centering is undefined " << centering)
774  return -1;
775  }
776  return 0;
777 }
778 
779 /** Given Cartesian mesh metadata extract whole_extent and bounds
780  * if bounds metadata is not already present then it is initialized
781  * from coordinate arrays. It's an error if whole_extent or coordinate
782  * arrays are not present. return zero if successful.
783  */
786  unsigned long *whole_extent, double *bounds);
787 
788 /// get the mesh's bounds from the coordinate axis arrays
792  double *bounds);
793 
794 /** Check that one Cartesian region covers the other coordinates must be in
795  * ascending order. assumes that both regions are specified in ascending order.
796  */
797 template <typename num_t>
799 int covers_ascending(const num_t *whole, const num_t *part)
800 {
801  if ((part[0] >= whole[0]) && (part[0] <= whole[1]) &&
802  (part[1] >= whole[0]) && (part[1] <= whole[1]) &&
803  (part[2] >= whole[2]) && (part[2] <= whole[3]) &&
804  (part[3] >= whole[2]) && (part[3] <= whole[3]) &&
805  (part[4] >= whole[4]) && (part[4] <= whole[5]) &&
806  (part[5] >= whole[4]) && (part[5] <= whole[5]))
807  return 1;
808  return 0;
809 }
810 
811 /** Check that one Cartesian region covers the other, taking into account the
812  * order of the coordinates. assumes that the regions are specified in the same
813  * orientation.
814  */
815 template <typename num_t>
817 int covers(const num_t *whole, const num_t *part)
818 {
819  bool x_ascend = whole[0] <= whole[1];
820  bool y_ascend = whole[2] <= whole[3];
821  bool z_ascend = whole[4] <= whole[5];
822  if (((x_ascend &&
823  (part[0] >= whole[0]) && (part[0] <= whole[1]) &&
824  (part[1] >= whole[0]) && (part[1] <= whole[1])) ||
825  (!x_ascend &&
826  (part[0] <= whole[0]) && (part[0] >= whole[1]) &&
827  (part[1] <= whole[0]) && (part[1] >= whole[1]))) &&
828  ((y_ascend &&
829  (part[2] >= whole[2]) && (part[2] <= whole[3]) &&
830  (part[3] >= whole[2]) && (part[3] <= whole[3])) ||
831  (!y_ascend &&
832  (part[2] <= whole[2]) && (part[2] >= whole[3]) &&
833  (part[3] <= whole[2]) && (part[3] >= whole[3]))) &&
834  ((z_ascend &&
835  (part[4] >= whole[4]) && (part[4] <= whole[5]) &&
836  (part[5] >= whole[4]) && (part[5] <= whole[5])) ||
837  (!z_ascend &&
838  (part[4] <= whole[4]) && (part[4] >= whole[5]) &&
839  (part[5] <= whole[4]) && (part[5] >= whole[5]))))
840  return 1;
841  return 0;
842 }
843 
844 /** check that two Cartesian regions have the same orientation ie they are
845  * either both specified in ascending or descending order.
846  */
847 template <typename num_t>
849 int same_orientation(const num_t *whole, const num_t *part)
850 {
851  if ((((whole[0] <= whole[1]) && (part[0] <= part[1])) ||
852  ((whole[0] >= whole[1]) && (part[0] >= part[1]))) &&
853  (((whole[2] <= whole[3]) && (part[2] <= part[3])) ||
854  ((whole[2] >= whole[3]) && (part[2] >= part[3]))) &&
855  (((whole[4] <= whole[5]) && (part[4] <= part[5])) ||
856  ((whole[4] >= whole[5]) && (part[4] >= part[5]))))
857  return 1;
858  return 0;
859 }
860 
861 /** where array dimensions specified by nx_max, ny_max, and nz_max are 1, and
862  * the extent would be out of bounds, set the extent to [0, 0]. If verbose is
863  * set, a warning is reported when the extent was clamped in one or more
864  * directions. The return is non zero if any direction was clamped and 0
865  * otherwise.
866  */
868 int clamp_dimensions_of_one(unsigned long nx_max, unsigned long ny_max,
869  unsigned long nz_max, unsigned long *extent, bool verbose);
870 
871 /** Return 0 if the passed extent does not exceed array dimensions specified in
872  * nx_max, ny_max, and nz_max. If verbose is set, an error is reported via
873  * TECA_ERROR when the extent would be out of bounds.
874  */
876 int validate_extent(unsigned long nx_max, unsigned long ny_max,
877  unsigned long nz_max, unsigned long *extent, bool verbose);
878 
879 
880 /// compares a set of arrays against a reference array
882 {
883 public:
884  /// set the array to which others will be compared to
885  void set_reference_array(const std::string &a_source,
886  const std::string &a_name, const std::string a_units,
887  const const_p_teca_variant_array &a_array);
888 
889  /// add an array to check against the reference
890  void append_array(const std::string &a_source,
891  const std::string &a_name, const std::string &a_units,
892  const const_p_teca_variant_array &a_array);
893 
894  /// error codes potentially returned from ::validate
895  enum {no_error = 0,
896  invalid_value = 1,
897  length_missmatch = 2,
898  type_missmatch = 3,
899  value_missmatch = 4,
900  unsupported_type = 5,
901  units_missmatch = 6
902  };
903 
904  /** Compare all the arrays in the collection against the reference returns
905  * 0 if all arrays in the collection are equal to the reference. When an
906  * array does not compare equal to the reference array a descritpion
907  * explaining why is returned in errorStr.
908  */
909  int validate(const std::string &a_descriptor,
910  double a_abs_tol, double a_rel_tol,
911  std::string &errorStr);
912 
913 private:
914  const_p_teca_variant_array m_reference;
915  std::string m_reference_source;
916  std::string m_reference_name;
917  std::string m_reference_units;
918 
919  std::vector<const_p_teca_variant_array> m_arrays;
920  std::vector<std::string> m_array_sources;
921  std::vector<std::string> m_array_names;
922  std::vector<std::string> m_array_units;
923 };
924 
925 /// Check that cooridnate arrays from different sources match a refrence array
926 /** Compares names, units, and values of coordinate axis arrays.
927  */
929 {
930 public:
932  m_absolute_tolerance(equal_tt<float>::absTol()),
933  m_relative_tolerance(equal_tt<float>::relTol())
934  {}
935 
936  /** Adds a time axis to validate. if provides_time is true then the axis
937  * becomes the reference to which all others are compared. returns 0 if
938  * necessary metadata was present and non zero if the necessary information
939  * was not present.
940  */
941  int add_time_axis(const std::string &source,
942  const teca_metadata &coords, const teca_metadata &atts,
943  bool provides_time);
944 
945  /** Adds anx-coordinate axis to validate. if provides_geometry is true then
946  * the axis becomes the reference to which all others are compared. returns
947  * 0 if necessary metadata was present and non zero if the necessary
948  * information was not present.
949  */
950  int add_x_coordinate_axis(const std::string &source,
951  const teca_metadata &coords, const teca_metadata &atts,
952  bool provides_geometry);
953 
954  /** Adds anx-coordinate axis to validate. if provides_geometry is true then
955  * the axis becomes the reference to which all others are compared. returns
956  * 0 if necessary metadata was present and non zero if the necessary
957  * information was not present.
958  */
959  int add_y_coordinate_axis(const std::string &source,
960  const teca_metadata &coords, const teca_metadata &atts,
961  bool provides_geometry);
962 
963  /** Adds anx-coordinate axis to validate. if provides_geometry is true then
964  * the axis becomes the reference to which all others are compared. returns
965  * 0 if necessary metadata was present and non zero if the necessary
966  * information was not present.
967  */
968  int add_z_coordinate_axis(const std::string &source,
969  const teca_metadata &coords, const teca_metadata &atts,
970  bool provides_geometry);
971 
972  /** runs the validation. returns 0 if all of the stored coordinate axes are
973  * equal to the reference axes. When an array does not compare equal to the
974  * reference array a descritpion explaining why is returned in errorStr.
975  */
976  int validate_spatial_coordinate_axes(std::string &errorStr);
977 
978  /** runs the validation. returns 0 if all of the stored time axes are
979  * equal to the reference axis. When an array does not compare equal to the
980  * reference array a descritpion explaining why is returned in errorStr.
981  */
982  int validate_time_axis(std::string &errorStr);
983 
984 private:
985  double m_absolute_tolerance;
986  double m_relative_tolerance;
991 };
992 
993 };
994 #endif
teca_coordinate_util::teca_validate_arrays
compares a set of arrays against a reference array
Definition: teca_coordinate_util.h:881
teca_variant_array.h
teca_coordinate_util::ascend_bracket
comparator implementing bracket for ascending input arrays
Definition: teca_coordinate_util.h:210
teca_variant_array_impl::get
void get(size_t i, U &val) const
Definition: teca_variant_array_impl.h:444
teca_metadata
A generic container for meta data in the form of name=value pairs.
Definition: teca_metadata.h:21
teca_coordinate_util::get_cartesian_mesh_bounds
TECA_EXPORT int get_cartesian_mesh_bounds(const const_p_teca_variant_array x, const const_p_teca_variant_array y, const const_p_teca_variant_array z, double *bounds)
get the mesh's bounds from the coordinate axis arrays
teca_coordinate_util::interpolate_t
struct TECA_EXPORT interpolate_t
A functor templated on order of accuracy for above Cartesian mesh interpolants.
Definition: teca_coordinate_util.h:677
teca_coordinate_util::time_step_of
TECA_EXPORT int time_step_of(const const_p_teca_variant_array &time, bool lower, bool clamp, const std::string &calendar, const std::string &units, const std::string &date, unsigned long &step)
teca_coordinate_util::bounds_to_extent
TECA_EXPORT int bounds_to_extent(const double *bounds, const const_p_teca_variant_array &x, const const_p_teca_variant_array &y, const const_p_teca_variant_array &z, unsigned long *extent)
teca_variant_array_impl
The concrete implementation of our type agnostic container for contiguous arrays.
Definition: teca_variant_array_impl.h:39
teca_coordinate_util
Codes dealing with operations on coordinate systems.
Definition: teca_coordinate_util.h:30
teca_coordinate_util::get_cartesian_mesh_extent
TECA_EXPORT int get_cartesian_mesh_extent(const teca_metadata &md, unsigned long *whole_extent, double *bounds)
teca_calcalcs::date
TECA_EXPORT int date(double val, int *year, int *month, int *day, int *hour, int *minute, double *second, const char *dataunits, const char *calendar_name)
teca_coordinate_util::validate_extent
TECA_EXPORT int validate_extent(unsigned long nx_max, unsigned long ny_max, unsigned long nz_max, unsigned long *extent, bool verbose)
teca_coordinate_util::equal
bool equal(T a, T b, T relTol=equal_tt< T >::relTol(), T absTol=equal_tt< T >::absTol(), typename std::enable_if< std::is_floating_point< T >::value >::type *=0)
Definition: teca_coordinate_util.h:70
teca_coordinate_util::get_table_offsets
TECA_EXPORT void get_table_offsets(const int_t *index, unsigned long n_rows, unsigned long &n_entities, std::vector< unsigned long > &counts, std::vector< unsigned long > &offsets, std::vector< unsigned long > &ids)
Definition: teca_coordinate_util.h:442
const_p_teca_variant_array
std::shared_ptr< const teca_variant_array > const_p_teca_variant_array
Definition: teca_variant_array.h:27
teca_coordinate_util::interpolate_linear
TECA_EXPORT int interpolate_linear(CT cx, CT cy, CT cz, const CT *p_x, const CT *p_y, const CT *p_z, const DT *p_data, unsigned long ihi, unsigned long jhi, unsigned long khi, unsigned long nx, unsigned long nxy, DT &val)
Definition: teca_coordinate_util.h:581
teca_coordinate_util::equal_tt
traits classes used to get default tolerances for comparing numbers of a given precision.
Definition: teca_coordinate_util.h:45
teca_coordinate_util::gt
Greater than predicate.
Definition: teca_coordinate_util.h:205
teca_coordinate_util::same_orientation
TECA_EXPORT int same_orientation(const num_t *whole, const num_t *part)
Definition: teca_coordinate_util.h:849
teca_coordinate_util::time_to_string
TECA_EXPORT int time_to_string(double val, const std::string &calendar, const std::string &units, const std::string &format, std::string &date)
teca_coordinate_util::interpolate_nearest
TECA_EXPORT int interpolate_nearest(CT cx, CT cy, CT cz, const CT *p_x, const CT *p_y, const CT *p_z, const DT *p_data, unsigned long ihi, unsigned long jhi, unsigned long khi, unsigned long nx, unsigned long nxy, DT &val)
Definition: teca_coordinate_util.h:492
teca_coordinate_util::equal_error
Definition: teca_coordinate_util.h:166
teca_coordinate_util::clamp_dimensions_of_one
TECA_EXPORT int clamp_dimensions_of_one(unsigned long nx_max, unsigned long ny_max, unsigned long nz_max, unsigned long *extent, bool verbose)
teca_coordinate_util::covers_ascending
TECA_EXPORT int covers_ascending(const num_t *whole, const num_t *part)
Definition: teca_coordinate_util.h:799
teca_coordinate_util::index_of
TECA_EXPORT int index_of(const data_t *data, unsigned long l, unsigned long r, data_t val, bool lower, unsigned long &id)
Definition: teca_coordinate_util.h:260
TEMPLATE_DISPATCH_FP
#define TEMPLATE_DISPATCH_FP(t, p, body)
Executes the code in body if p is a t<nt> where nt is a floating point type.
Definition: teca_variant_array_impl.h:171
teca_coordinate_util::descend_bracket
comparator implementing bracket for descending input arrays
Definition: teca_coordinate_util.h:232
max_prec
#define max_prec(T)
For printing data as ASCII with the maximum supported numerical precision.
Definition: teca_coordinate_util.h:26
teca_coordinate_util::teca_coordinate_axis_validator
Check that cooridnate arrays from different sources match a refrence array.
Definition: teca_coordinate_util.h:928
teca_coordinate_util::covers
TECA_EXPORT int covers(const num_t *whole, const num_t *part)
Definition: teca_coordinate_util.h:817
teca_coordinate_util::geq
Greater than or equal to predicate.
Definition: teca_coordinate_util.h:195
teca_coordinate_util::convert_cell_extent
TECA_EXPORT int convert_cell_extent(num_t *extent, int centering)
convert from a cell extent to a face, edge or point centered extent
Definition: teca_coordinate_util.h:734
teca_error::TECA_EXPORT
p_teca_error_handler error_handler TECA_EXPORT
The global error handler instance.
teca_coordinate_util::lt
Less than predicate.
Definition: teca_coordinate_util.h:200
TECA_ERROR
#define TECA_ERROR(_msg)
Constructs an error message and sends it to the stderr stream.
Definition: teca_common.h:146
teca_coordinate_util::validate_centering
TECA_EXPORT int validate_centering(int centering)
return 0 if the centering is one of the values defined in teca_array_attributes
teca_coordinate_util::leq
Less than or equal to predicate.
Definition: teca_coordinate_util.h:190
teca_variant_array_impl.h