nidas  v1.2-1520
Wind3D.h
Go to the documentation of this file.
1 // -*- mode: C++; indent-tabs-mode: nil; c-basic-offset: 4; tab-width: 4; -*-
2 // vim: set shiftwidth=4 softtabstop=4 expandtab:
3 /*
4  ********************************************************************
5  ** NIDAS: NCAR In-situ Data Acquistion Software
6  **
7  ** 2006, Copyright University Corporation for Atmospheric Research
8  **
9  ** This program is free software; you can redistribute it and/or modify
10  ** it under the terms of the GNU General Public License as published by
11  ** the Free Software Foundation; either version 2 of the License, or
12  ** (at your option) any later version.
13  **
14  ** This program is distributed in the hope that it will be useful,
15  ** but WITHOUT ANY WARRANTY; without even the implied warranty of
16  ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  ** GNU General Public License for more details.
18  **
19  ** The LICENSE.txt file accompanying this software contains
20  ** a copy of the GNU General Public License. If it is not found,
21  ** write to the Free Software Foundation, Inc.,
22  ** 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
23  **
24  ********************************************************************
25 */
26 
27 #ifndef NIDAS_DNYLD_ISFF_WIND3D_H
28 #define NIDAS_DNYLD_ISFF_WIND3D_H
29 
32 #include <nidas/Config.h>
33 #include "WindOrienter.h"
34 #include "WindTilter.h"
35 #include "WindRotator.h"
36 
37 #ifdef HAVE_LIBGSL
38 #include <gsl/gsl_linalg.h>
39 #endif
40 
41 namespace nidas {
42 
43 namespace dynld { namespace isff {
44 
50 {
51 public:
52 
53  Wind3D();
54 
55  ~Wind3D();
56 
62  bool process(const nidas::core::Sample* samp,
63  std::list<const nidas::core::Sample*>& results) throw();
64 
65  void setBias(int i, double val);
66 
67  double getBias(int i) const
68  {
69  return _bias[i];
70  }
71 
72  double getVazimuth() const
73  {
74  return _rotator.getAngleDegrees();
75  }
76 
87  void setVazimuth(double val)
88  {
90  }
91 
92  double getLeanDegrees() const
93  {
94  return _tilter.getLeanDegrees();
95  }
96  void setLeanDegrees(double val)
97  {
99  }
100 
101  double getLeanAzimuthDegrees() const
102  {
104  }
105  void setLeanAzimuthDegrees(double val)
106  {
108  }
109 
110  void setDespike(bool val)
111  {
112  _despike = val;
113  }
114 
115  bool getDespike() const
116  {
117  return _despike;
118  }
119 
120  void setMetek(int ismetek) {
121  _metek = (ismetek == 1);
122  }
123 
124  void setOutlierProbability(double val)
125  {
126  for (int i = 0; i < 4; i++)
128  }
129 
130  double getOutlierProbability() const
131  {
132  return _despiker[0].getOutlierProbability();
133  }
134 
135  void setDiscLevelMultiplier(double val)
136  {
137  for (int i = 0; i < 4; i++)
139  }
140 
141  double getDiscLevelMultiplier() const
142  {
143  return _despiker[0].getDiscLevelMultiplier();
144  }
145 
146  double getDiscLevel() const
147  {
148  return _despiker[0].getDiscLevel();
149  }
150 
151  void setTcOffset(double val)
152  {
153  _tcOffset = val;
154  }
155 
156  double getTcOffset() const
157  {
158  return _tcOffset;
159  }
160 
161  void setTcSlope(double val)
162  {
163  _tcSlope = val;
164  }
165 
166  double getTcSlope() const
167  {
168  return _tcSlope;
169  }
170 
174  void setDoHorizontalRotation(bool val)
175  {
176  _horizontalRotation = val;
177  }
178 
182  void setDoTiltCorrection(bool val)
183  {
184  _tiltCorrection = val;
185  }
186 
187  void despike(nidas::core::dsm_time_t tt,float* uvwt,int n, bool* spikeOrMissing)
188  throw();
189 
200  void offsetsTiltAndRotate(nidas::core::dsm_time_t tt, float* uvwt) throw();
201 
205  void applyOrientation(nidas::core::dsm_time_t tt, float* uvwt) throw();
206 
212 
217  void validate()
218  throw(nidas::util::InvalidParameterException);
219 
224  void validateSscanfs() throw(nidas::util::InvalidParameterException);
225 
233  virtual void parseParameters() throw(nidas::util::InvalidParameterException);
234 
238  virtual void checkSampleTags() throw(nidas::util::InvalidParameterException);
239 
240 #ifdef HAVE_LIBGSL
241 
246  virtual void getTransducerRotation(nidas::core::dsm_time_t tt) throw();
247 
248  virtual void transducerShadowCorrection(nidas::core::dsm_time_t, float *) throw();
249 #endif
250 
251 protected:
252 
256 
257  static const int DATA_GAP_USEC = 60000000;
258 
260 
261  double _bias[3];
262 
264 
265  bool _despike;
266 
267  bool _metek;
268 
270 
272 
274 
276 
277  double _tcOffset;
278 
279  double _tcSlope;
280 
285 
290 
295 
301 
307 
313 
319 
320  unsigned int _noutVals;
321 
326  unsigned int _numParsed;
327 
332 
333 #ifdef HAVE_LIBGSL
334 
339  nidas::core::CalFile* _atCalFile;
340 
344  double _atMatrix[3][3];
345 
346 #define COMPUTE_ABC2UVW_INVERSE
347 #ifdef COMPUTE_ABC2UVW_INVERSE
348  double _atInverse[3][3];
349 #else
350  gsl_vector* _atVectorGSL1;
351  gsl_vector* _atVectorGSL2;
352 #endif
353 
354  gsl_matrix* _atMatrixGSL;
355 
356  gsl_permutation* _atPermutationGSL;
357 #endif
358 
365 
366 private:
367 
368  // no copying
369  Wind3D(const Wind3D& x);
370 
371  // no assignment
372  Wind3D& operator=(const Wind3D& x);
373 
374 };
375 
376 }}} // namespace nidas namespace dynld namespace isff
377 
378 #endif
double getTcSlope() const
Definition: Wind3D.h:166
A class for reading ASCII files containing a time series of calibration data.
Definition: CalFile.h:164
double getOutlierProbability() const
Definition: Wind3D.h:130
void setTcOffset(double val)
Definition: Wind3D.h:151
A class for performing the common processes on wind data from a 3D sonic anemometer.
Definition: Wind3D.h:49
void despike(nidas::core::dsm_time_t tt, float *uvwt, int n, bool *spikeOrMissing)
Definition: Wind3D.cc:97
nidas::dynld::isff::WindTilter WindTilter
Definition: Wind3D.h:255
bool _horizontalRotation
Should horizontal rotation of U,V be performed?
Definition: Wind3D.h:284
virtual void checkSampleTags()
Check the SampleTags that are defined for this sensor.
Definition: Wind3D.cc:357
unsigned int dsm_sample_id_t
Definition: Sample.h:63
unsigned int _noutVals
Definition: Wind3D.h:320
Adaptive forecaster for despiking of time-series data.
Definition: AdaptiveDespiker.h:41
nidas::core::CalFile * _oaCalFile
CalFile containing wind offsets and rotation angles.
Definition: Wind3D.h:331
Support for a sensor that is sending packets on a TCP socket, a UDP socket, a Bluetooth RF Comm socke...
Definition: SerialSensor.h:64
static const int DATA_GAP_USEC
Definition: Wind3D.h:257
bool _metek
Definition: Wind3D.h:267
void validateSscanfs()
Warn user if number of scanf fields does not match number expected from variables in sample...
Definition: Wind3D.cc:400
void setLeanDegrees(double val)
Definition: WindTilter.h:129
long long dsm_time_t
Posix time in microseconds, the number of non-leap microseconds since 1970 Jan 1 00:00 UTC...
Definition: Sample.h:61
void applyOrientation(nidas::core::dsm_time_t tt, float *uvwt)
Apply orientation changes to the wind components.
Definition: Wind3D.cc:205
void setLeanAzimuthDegrees(double val)
Definition: Wind3D.h:105
float getDiscLevelMultiplier() const
Definition: AdaptiveDespiker.h:70
Wind3D & operator=(const Wind3D &x)
void setVazimuth(double val)
Wind vectors in geographic coordinates are expressed by U, the component of the wind blowing toward t...
Definition: Wind3D.h:87
void validate()
Validate the configuration of this sensor.
Definition: Wind3D.cc:225
double getLeanAzimuthDegrees() const
Definition: WindTilter.h:135
WindRotator _rotator
Definition: Wind3D.h:271
bool _tiltCorrection
Should 3D tilt correction be applied?
Definition: Wind3D.h:289
void setDoTiltCorrection(bool val)
Should 3D tilt corrections be applied?
Definition: Wind3D.h:182
WindTilter is used to apply a 3d rotation matrix to a wind vector.
Definition: WindTilter.h:119
void setTcSlope(double val)
Definition: Wind3D.h:161
void setLeanDegrees(double val)
Definition: Wind3D.h:96
double getTcOffset() const
Definition: Wind3D.h:156
virtual void parseParameters()
Parse the list of nidas::core::Parameter that are associated with this sensor.
Definition: Wind3D.cc:234
void setOutlierProbability(double val)
Definition: Wind3D.h:124
double getVazimuth() const
Definition: Wind3D.h:72
void offsetsTiltAndRotate(nidas::core::dsm_time_t tt, float *uvwt)
Do standard bias removal, tilt correction and horizontal rotation of 3d sonic anemometer data...
Definition: Wind3D.cc:211
unsigned int _numParsed
Number of variables that are parsed from input, i.e.
Definition: Wind3D.h:326
double getLeanAzimuthDegrees() const
Definition: Wind3D.h:101
void setMetek(int ismetek)
Definition: Wind3D.h:120
double getAngleDegrees() const
Definition: WindRotator.cc:14
nidas::core::AdaptiveDespiker _despiker[4]
Definition: Wind3D.h:269
Wind3D()
Definition: Wind3D.cc:51
void setAngleDegrees(double val)
Definition: WindRotator.cc:19
nidas::dynld::isff::WindRotator WindRotator
Definition: Wind3D.h:254
WindOrienter _orienter
Definition: Wind3D.h:275
float getOutlierProbability() const
Definition: AdaptiveDespiker.h:59
double _tcSlope
Definition: Wind3D.h:279
Interface to a data sample.
Definition: Sample.h:189
int _ldiagIndex
If user requests &quot;ldiag&quot;, its index in the output sample.
Definition: Wind3D.h:306
void setDoHorizontalRotation(bool val)
Should 2D horizontal rotations of U,V be applied?
Definition: Wind3D.h:174
double _shadowFactor
Transducer shadow (aka flow distortion) correction factor.
Definition: Wind3D.h:364
double _bias[3]
Definition: Wind3D.h:261
bool _despike
Definition: Wind3D.h:265
double _tcOffset
Definition: Wind3D.h:277
int _diagIndex
If user requests &quot;diag&quot; or &quot;status&quot;, its index in the output sample.
Definition: Wind3D.h:300
double getLeanDegrees() const
Definition: Wind3D.h:92
bool getDespike() const
Definition: Wind3D.h:115
nidas::core::dsm_time_t _ttlast[4]
Definition: Wind3D.h:259
double getLeanDegrees() const
Definition: WindTilter.h:124
void setDespike(bool val)
Definition: Wind3D.h:110
void readOffsetsAnglesCalFile(nidas::core::dsm_time_t tt)
Update the settings from the offsets and angles calibration file, if any.
Definition: Wind3D.cc:136
bool process(const nidas::core::Sample *samp, std::list< const nidas::core::Sample * > &results)
Basic process method for sonic anemometer wind plus temperature: u,v,w,tc parsed from an ASCII sample...
Definition: Wind3D.cc:543
bool _allBiasesNaN
Definition: Wind3D.h:263
nidas::core::dsm_sample_id_t _sampleId
Id of output sample.
Definition: Wind3D.h:294
A class for rotating winds according to different orientations of the wind sensor.
Definition: WindOrienter.h:44
int _dirIndex
If user requests wind direction, variable name &quot;dir&quot;, its index in the output sample.
Definition: Wind3D.h:318
~Wind3D()
Definition: Wind3D.cc:86
nidas::dynld::isff::WindOrienter WindOrienter
Definition: Wind3D.h:253
void setLeanAzimuthDegrees(double val)
Definition: WindTilter.h:140
void setDiscLevelMultiplier(double val)
Definition: Wind3D.h:135
float getDiscLevel() const
If a value is more than discrLevel * sigma away from the mean, then it is considered a spike...
Definition: AdaptiveDespiker.h:77
void setBias(int i, double val)
Definition: Wind3D.cc:124
int _spdIndex
If user requests wind speed, variable name &quot;spd&quot;, its index in the output sample. ...
Definition: Wind3D.h:312
double getBias(int i) const
Definition: Wind3D.h:67
double getDiscLevelMultiplier() const
Definition: Wind3D.h:141
WindTilter _tilter
Definition: Wind3D.h:273
Rotate a (U,V) 2D wind vector by an angle.
Definition: WindRotator.h:39
Definition: InvalidParameterException.h:35
double getDiscLevel() const
Definition: Wind3D.h:146