nidas v1.2.3
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
41namespace nidas {
42
43namespace dynld { namespace isff {
44
50{
51public:
52
53 Wind3D();
54
55 ~Wind3D();
56
62 bool process(const nidas::core::Sample* samp,
63 std::list<const nidas::core::Sample*>& results);
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 {
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
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
131 {
133 }
134
135 void setDiscLevelMultiplier(double val)
136 {
137 for (int i = 0; i < 4; i++)
139 }
140
142 {
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
175 {
177 }
178
182 void setDoTiltCorrection(bool val)
183 {
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
223 void validateSscanfs();
224
232 virtual void parseParameters();
233
237 virtual void checkSampleTags();
238
239#ifdef HAVE_LIBGSL
246
247 virtual void transducerShadowCorrection(nidas::core::dsm_time_t, float *);
248#endif
249
250protected:
251
255
256 static const int DATA_GAP_USEC = 60000000;
257
259
260 double _bias[3];
261
263
265
266 bool _metek;
267
269
271
273
275
276 double _tcOffset;
277
278 double _tcSlope;
279
284
289
294
300
306
312
318
319 unsigned int _noutVals;
320
325 unsigned int _numParsed;
326
331
332#ifdef HAVE_LIBGSL
339
343 double _atMatrix[3][3];
344
345#define COMPUTE_ABC2UVW_INVERSE
346#ifdef COMPUTE_ABC2UVW_INVERSE
347 double _atInverse[3][3];
348#else
349 gsl_vector* _atVectorGSL1;
350 gsl_vector* _atVectorGSL2;
351#endif
352
353 gsl_matrix* _atMatrixGSL;
354
355 gsl_permutation* _atPermutationGSL;
356#endif
357
364
365private:
366
367 // no copying
368 Wind3D(const Wind3D& x);
369
370 // no assignment
372
373};
374
375}}} // namespace nidas namespace dynld namespace isff
376
377#endif
Adaptive forecaster for despiking of time-series data.
Definition AdaptiveDespiker.h:42
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
float getDiscLevelMultiplier() const
Definition AdaptiveDespiker.h:70
float getOutlierProbability() const
Definition AdaptiveDespiker.h:59
A class for reading ASCII files containing a time series of calibration data.
Definition CalFile.h:166
Interface to a data sample.
Definition Sample.h:190
Support for a sensor that is sending packets on a TCP socket, a UDP socket, a Bluetooth RF Comm socke...
Definition SerialSensor.h:65
A class for performing the common processes on wind data from a 3D sonic anemometer.
Definition Wind3D.h:50
void despike(nidas::core::dsm_time_t tt, float *uvwt, int n, bool *spikeOrMissing)
Definition Wind3D.cc:97
void setDiscLevelMultiplier(double val)
Definition Wind3D.h:135
unsigned int _noutVals
Definition Wind3D.h:319
void applyOrientation(nidas::core::dsm_time_t tt, float *uvwt)
Apply orientation changes to the wind components.
Definition Wind3D.cc:205
nidas::core::dsm_sample_id_t _sampleId
Id of output sample.
Definition Wind3D.h:293
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
double _bias[3]
Definition Wind3D.h:260
nidas::core::dsm_time_t _ttlast[4]
Definition Wind3D.h:258
double getTcOffset() const
Definition Wind3D.h:156
Wind3D()
Definition Wind3D.cc:51
void setTcOffset(double val)
Definition Wind3D.h:151
double getLeanDegrees() const
Definition Wind3D.h:92
void setBias(int i, double val)
Definition Wind3D.cc:124
double getLeanAzimuthDegrees() const
Definition Wind3D.h:101
void setDoTiltCorrection(bool val)
Should 3D tilt corrections be applied?
Definition Wind3D.h:182
bool _horizontalRotation
Should horizontal rotation of U,V be performed?
Definition Wind3D.h:283
int _spdIndex
If user requests wind speed, variable name "spd", its index in the output sample.
Definition Wind3D.h:311
int _ldiagIndex
If user requests "ldiag", its index in the output sample.
Definition Wind3D.h:305
double getVazimuth() const
Definition Wind3D.h:72
nidas::core::CalFile * _oaCalFile
CalFile containing wind offsets and rotation angles.
Definition Wind3D.h:330
double getDiscLevel() const
Definition Wind3D.h:146
void setLeanDegrees(double val)
Definition Wind3D.h:96
void setLeanAzimuthDegrees(double val)
Definition Wind3D.h:105
double _shadowFactor
Transducer shadow (aka flow distortion) correction factor.
Definition Wind3D.h:363
nidas::dynld::isff::WindRotator WindRotator
Definition Wind3D.h:253
bool getDespike() const
Definition Wind3D.h:115
int _diagIndex
If user requests "diag" or "status", its index in the output sample.
Definition Wind3D.h:299
nidas::dynld::isff::WindOrienter WindOrienter
Definition Wind3D.h:252
double getBias(int i) const
Definition Wind3D.h:67
double getTcSlope() const
Definition Wind3D.h:166
Wind3D(const Wind3D &x)
void setTcSlope(double val)
Definition Wind3D.h:161
void validate()
Validate the configuration of this sensor.
Definition Wind3D.cc:225
bool _metek
Definition Wind3D.h:266
double _tcOffset
Definition Wind3D.h:276
double _tcSlope
Definition Wind3D.h:278
double getOutlierProbability() const
Definition Wind3D.h:130
nidas::core::AdaptiveDespiker _despiker[4]
Definition Wind3D.h:268
int _dirIndex
If user requests wind direction, variable name "dir", its index in the output sample.
Definition Wind3D.h:317
bool _despike
Definition Wind3D.h:264
void setDoHorizontalRotation(bool val)
Should 2D horizontal rotations of U,V be applied?
Definition Wind3D.h:174
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
Wind3D & operator=(const Wind3D &x)
static const int DATA_GAP_USEC
Definition Wind3D.h:256
virtual void checkSampleTags()
Check the SampleTags that are defined for this sensor.
Definition Wind3D.cc:356
bool _tiltCorrection
Should 3D tilt correction be applied?
Definition Wind3D.h:288
double getDiscLevelMultiplier() const
Definition Wind3D.h:141
void setDespike(bool val)
Definition Wind3D.h:110
~Wind3D()
Definition Wind3D.cc:86
void setOutlierProbability(double val)
Definition Wind3D.h:124
bool _allBiasesNaN
Definition Wind3D.h:262
virtual void parseParameters()
Parse the list of nidas::core::Parameter that are associated with this sensor.
Definition Wind3D.cc:234
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:540
unsigned int _numParsed
Number of variables that are parsed from input, i.e.
Definition Wind3D.h:325
WindTilter _tilter
Definition Wind3D.h:272
void setMetek(int ismetek)
Definition Wind3D.h:120
WindRotator _rotator
Definition Wind3D.h:270
void readOffsetsAnglesCalFile(nidas::core::dsm_time_t tt)
Update the settings from the offsets and angles calibration file, if any.
Definition Wind3D.cc:136
void validateSscanfs()
Warn user if number of scanf fields does not match number expected from variables in sample.
Definition Wind3D.cc:397
WindOrienter _orienter
Definition Wind3D.h:274
nidas::dynld::isff::WindTilter WindTilter
Definition Wind3D.h:254
A class for rotating winds according to different orientations of the wind sensor.
Definition WindOrienter.h:45
Rotate a (U,V) 2D wind vector by an angle.
Definition WindRotator.h:39
void setAngleDegrees(double val)
Definition WindRotator.cc:19
double getAngleDegrees() const
Definition WindRotator.cc:14
WindTilter is used to apply a 3d rotation matrix to a wind vector.
Definition WindTilter.h:119
double getLeanAzimuthDegrees() const
Definition WindTilter.h:135
void setLeanDegrees(double val)
Definition WindTilter.h:129
void setLeanAzimuthDegrees(double val)
Definition WindTilter.h:140
double getLeanDegrees() const
Definition WindTilter.h:124
Sample * getSample(sampleType type, unsigned int len)
A convienence method for getting a sample of an enumerated type from a pool.
Definition Sample.cc:70
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:62
unsigned int dsm_sample_id_t
Definition Sample.h:64
Root namespace for the NCAR In-Situ Data Acquisition Software.
Definition A2DConverter.h:31