MeteoIODoc 20250312.660e6d76
InterpolARIMA.h
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-3.0-or-later
2/***********************************************************************************/
3/* Copyright 2013 WSL Institute for Snow and Avalanche Research SLF-DAVOS */
4/***********************************************************************************/
5/* This file is part of MeteoIO.
6 MeteoIO is free software: you can redistribute it and/or modify
7 it under the terms of the GNU Lesser General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 MeteoIO is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU Lesser General Public License for more details.
15
16 You should have received a copy of the GNU Lesser General Public License
17 along with MeteoIO. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef INTERPOLARIMA_H
21#define INTERPOLARIMA_H
22
23/* For peace of mind, disable strict ISO c++ compliance as it would emit warnings
24 * because of ctsa's variable lenght arrays*/
25#ifdef __GNUC__
26 #pragma GCC diagnostic push
27 #pragma GCC diagnostic ignored "-Wpedantic"
28 #include <meteoio/thirdParty/ctsa.h>
29 #pragma GCC diagnostic pop
30#elif defined __clang__
31 #pragma clang diagnostic push
32 #pragma clang diagnostic ignored "-Wpedantic"
33 #include <meteoio/thirdParty/ctsa.h>
34 #pragma clang diagnostic pop
35#else
36 #include <meteoio/thirdParty/ctsa.h>
37#endif
38
39#include <iostream>
40#include <map>
42#include <string>
43#include <vector>
44
45namespace mio {
46
47using namespace ARIMAutils;
105 public:
107 InterpolARIMA(const std::vector<double>& data_in, const size_t& gap_loc, const size_t& N_gap, const int& s = 0);
108 InterpolARIMA(const std::vector<double>& data_in, const size_t& gap_loc, const size_t& N_gap, const std::vector<double>& xreg_vec, const int& s = 0);
109 InterpolARIMA(const std::vector<double>& data_in, const size_t& gap_loc, const size_t& n_predictions, const std::string& direction = "forward", const int& s = 0);
110
111 // Setters
112 void setAutoArimaMetaData(int max_p_param = 8, int max_d_param = 3, int max_q = 8, int start_p = 2, int start_q = 2, int max_P = 2,
113 int max_D = 1, int max_Q = 2, int start_P = 1, int start_Q = 1, bool seasonal = true,
114 bool stationary = false);
115 void setOptMetaData(ObjectiveFunction method = CSS_MLE, OptimizationMethod opt_method = BFGS, bool stepwise = true,
116 bool approximation = false, int num_models = 94);
117 void setVerbose(bool verbose = false);
119 void setManualARIMA(int p, int d, int q, int P, int D, int Q, bool fill_backward);
120
121 // Interpolation methods
122 void fillGap();
123 void fillGapManual();
124
125 void interpolate();
126 std::vector<double> predict(size_t n_steps = 0);
127 // only do denormalization ARIMA predict! otherwise it will denorm twice
128 std::vector<double> ARIMApredict(size_t n_steps);
129
130 // Getters
131 std::vector<double> getData() { return norm.denormalize(data); }
132 std::vector<double> getForwardData() { return norm.denormalize(data_forward); }
133 std::vector<double> getBackwardData() { return norm.denormalize(data_backward); }
134 std::vector<double> getInterpolatedData();
135
136 // Copy constructor
138 : norm(other.norm), data(other.data), gap_loc(other.gap_loc), N_gap(other.N_gap), time(other.time), pred_forward(other.pred_forward),
139 pred_backward(other.pred_backward), xreg_vec_f(other.xreg_vec_f), xreg_vec_b(other.xreg_vec_b),
140 data_forward(other.data_forward), data_backward(other.data_backward), new_xreg_vec_f(other.new_xreg_vec_f),
141 new_xreg_vec_b(other.new_xreg_vec_b), xreg_f((xreg_vec_f.empty()) ? nullptr : &xreg_vec_f[0]),
142 xreg_b((xreg_vec_b.empty()) ? nullptr : &xreg_vec_b[0]), new_xreg_f((new_xreg_vec_f.empty()) ? nullptr : &new_xreg_vec_f[0]),
143 new_xreg_b((new_xreg_vec_b.empty()) ? nullptr : &new_xreg_vec_b[0]), amse_forward(other.amse_forward),
144 amse_backward(other.amse_backward), N_data_forward(other.N_data_forward), N_data_backward(other.N_data_backward),
145 max_p(other.max_p), max_d(other.max_d), max_q(other.max_q), start_p(other.start_p), start_q(other.start_q),
146 max_P(other.max_P), max_D(other.max_D), max_Q(other.max_Q), start_P(other.start_P), start_Q(other.start_Q), r(other.r),
147 s(other.s), method(other.method), opt_method(other.opt_method), stepwise(other.stepwise), approximation(other.approximation),
148 num_models(other.num_models), seasonal(other.seasonal), stationary(other.stationary),
149 auto_arima_forward(auto_arima_copy(other.auto_arima_forward)), auto_arima_backward(auto_arima_copy(other.auto_arima_backward)), sarima_forward(other.sarima_forward) {
150 }
151
152 // Copy assignment operator
154 operator=(const InterpolARIMA &other) {
155 // protect against invalid self-assignment
156 if (this != &other) {
157 auto_arima_forward = auto_arima_copy(other.auto_arima_forward);
158 auto_arima_backward = auto_arima_copy(other.auto_arima_backward);
159
160 // 3: copy all the other fields from the other object
161 gap_loc = other.gap_loc;
162 N_gap = other.N_gap;
163 time = other.time;
164 pred_forward = other.pred_forward;
165 pred_backward = other.pred_backward;
166 norm = other.norm;
167 data = other.data;
168 xreg_vec_f = other.xreg_vec_f;
169 xreg_vec_b = other.xreg_vec_b;
170 data_forward = other.data_forward;
171 data_backward = other.data_backward;
172 new_xreg_vec_f = other.new_xreg_vec_f;
173 new_xreg_vec_b = other.new_xreg_vec_b;
174 N_data_forward = other.N_data_forward;
175 N_data_backward = other.N_data_backward;
176 max_p = other.max_p;
177 max_d = other.max_d;
178 max_q = other.max_q;
179 start_p = other.start_p;
180 start_q = other.start_q;
181 max_P = other.max_P;
182 max_D = other.max_D;
183 max_Q = other.max_Q;
184 start_P = other.start_P;
185 start_Q = other.start_Q;
186 r = other.r;
187 s = other.s;
188 method = other.method;
189 opt_method = other.opt_method;
190 stepwise = other.stepwise;
191 approximation = other.approximation;
192 num_models = other.num_models;
193 seasonal = other.seasonal;
194 stationary = other.stationary;
195
196 // 4: handle the pointers to the vectors
197 xreg_f = (xreg_vec_f.empty()) ? nullptr : &xreg_vec_f[0];
198 xreg_b = (xreg_vec_b.empty()) ? nullptr : &xreg_vec_b[0];
199 new_xreg_f = (new_xreg_vec_f.empty()) ? nullptr : &new_xreg_vec_f[0];
200 new_xreg_b = (new_xreg_vec_b.empty()) ? nullptr : &new_xreg_vec_b[0];
201 }
202 // by convention, always return *this
203 return *this;
204 }
205
207 auto_arima_free(auto_arima_forward);
208 auto_arima_free(auto_arima_backward);
209 delete xreg_f;
210 delete xreg_b;
211 delete new_xreg_f;
212 delete new_xreg_b;
213 }
214
215 // info
216 std::string toString();
217 std::string autoArimaInfo(const auto_arima_object& obj);
218
219
220private:
221 // Interpolation variables
222 Normalization norm;
223 std::vector<double> data;
224 size_t gap_loc;
225 size_t N_gap;
226 std::vector<double> time;
227 std::vector<double> pred_forward, pred_backward;
228
229 // Auto Arima variables
230 // const doesnt work with c
231 std::vector<double> xreg_vec_f, xreg_vec_b, data_forward, data_backward, new_xreg_vec_f, new_xreg_vec_b;
232 double *xreg_f;
233 double *xreg_b;
234 double *new_xreg_f;
235 double *new_xreg_b;
236 std::vector<double> amse_forward, amse_backward;
237 size_t N_data_forward, N_data_backward;
238 int max_p = 8, max_d = 3, max_q = 8;
239 int start_p = 2, start_q = 2;
240 int max_P = 2, max_D = 1, max_Q = 2;
241 int start_P = 1, start_Q = 1;
242 int r = 0, s = 0;
243 ObjectiveFunction method = CSS_MLE; OptimizationMethod opt_method = BFGS;
244 bool stepwise = true, approximation = true;
245 int num_models = 94;
246 bool seasonal = true, stationary = false;
247
248 bool consistencyCheck();
249 auto_arima_object initAutoArima(const size_t& N_data);
250
251 // last to be initialized
252public:
253 auto_arima_object auto_arima_forward;
254 auto_arima_object auto_arima_backward;
255
256private:
257 // (S)ARIMA variables
258 bool set_manual = false;
259 bool fill_backward_manual = false;
260
261public:
262 sarima_object sarima_forward;
263};
264
265} // namespace mio
266
267#endif // INTERPOLARIMA_H
Definition: ARIMAutils.h:55
std::vector< double > denormalize(const std::vector< double > &data)
Definition: ARIMAutils.cc:80
Mode
Definition: ARIMAutils.h:57
This class is used for interpolating or predicting missing data in a time series using the Auto ARIMA...
Definition: InterpolARIMA.h:104
std::string toString()
Definition: InterpolARIMA.cc:129
sarima_object sarima_forward
Definition: InterpolARIMA.h:262
InterpolARIMA & operator=(const InterpolARIMA &other)
Definition: InterpolARIMA.h:154
InterpolARIMA()
Definition: InterpolARIMA.cc:34
std::vector< double > getBackwardData()
Definition: InterpolARIMA.h:133
std::vector< double > getForwardData()
Definition: InterpolARIMA.h:132
void fillGap()
Definition: InterpolARIMA.cc:420
~InterpolARIMA()
Definition: InterpolARIMA.h:206
void setAutoArimaMetaData(int max_p_param=8, int max_d_param=3, int max_q=8, int start_p=2, int start_q=2, int max_P=2, int max_D=1, int max_Q=2, int start_P=1, int start_Q=1, bool seasonal=true, bool stationary=false)
Definition: InterpolARIMA.cc:312
InterpolARIMA(const InterpolARIMA &other)
Definition: InterpolARIMA.h:137
auto_arima_object auto_arima_forward
Definition: InterpolARIMA.h:253
std::vector< double > ARIMApredict(size_t n_steps)
Definition: InterpolARIMA.cc:393
std::string autoArimaInfo(const auto_arima_object &obj)
Definition: InterpolARIMA.cc:175
std::vector< double > getData()
Definition: InterpolARIMA.h:131
void setManualARIMA(int p, int d, int q, int P, int D, int Q, bool fill_backward)
Definition: InterpolARIMA.cc:375
void interpolate()
Definition: InterpolARIMA.cc:535
void fillGapManual()
Definition: InterpolARIMA.cc:463
void setNormalizationMode(Normalization::Mode mode)
Definition: InterpolARIMA.cc:298
std::vector< double > getInterpolatedData()
Definition: InterpolARIMA.cc:386
void setOptMetaData(ObjectiveFunction method=CSS_MLE, OptimizationMethod opt_method=BFGS, bool stepwise=true, bool approximation=false, int num_models=94)
Definition: InterpolARIMA.cc:356
std::vector< double > predict(size_t n_steps=0)
Definition: InterpolARIMA.cc:560
void setVerbose(bool verbose=false)
Definition: InterpolARIMA.cc:370
auto_arima_object auto_arima_backward
Definition: InterpolARIMA.h:254
OptimizationMethod
Definition: ARIMAutils.h:41
@ BFGS
Definition: ARIMAutils.h:47
ObjectiveFunction
Definition: ARIMAutils.h:34
@ CSS_MLE
Definition: ARIMAutils.h:35
Definition: Config.cc:31