MeteoIODoc 20241221.207bde49
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(std::vector<double> data_in, size_t gap_loc, size_t N_gap, int s = 0);
108 InterpolARIMA(std::vector<double> data_in, size_t gap_loc, size_t N_gap, std::vector<double> xreg_vec, int s = 0);
109 InterpolARIMA(std::vector<double> data_in, size_t gap_loc, size_t n_predictions, std::string direction = "forward", 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 std::vector<double> simulate(int n_steps, int seed = 0);
123 void fillGap();
124 void fillGapManual();
125
126 void interpolate();
127 std::vector<double> predict(size_t n_steps = 0);
128 // only do denormalization ARIMA predict! otherwise it will denorm twice
129 std::vector<double> ARIMApredict(size_t n_steps);
130
131 // Getters
132 std::vector<double> getData() { return norm.denormalize(data); }
133 std::vector<double> getForwardData() { return norm.denormalize(data_forward); }
134 std::vector<double> getBackwardData() { return norm.denormalize(data_backward); }
135 std::vector<double> getInterpolatedData();
136
137 // Copy constructor
139 : norm(other.norm), data(other.data), gap_loc(other.gap_loc), N_gap(other.N_gap), time(other.time), pred_forward(other.pred_forward),
140 pred_backward(other.pred_backward), xreg_vec_f(other.xreg_vec_f), xreg_vec_b(other.xreg_vec_b),
141 data_forward(other.data_forward), data_backward(other.data_backward), new_xreg_vec_f(other.new_xreg_vec_f),
142 new_xreg_vec_b(other.new_xreg_vec_b), xreg_f((xreg_vec_f.empty()) ? nullptr : &xreg_vec_f[0]),
143 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]),
144 new_xreg_b((new_xreg_vec_b.empty()) ? nullptr : &new_xreg_vec_b[0]), amse_forward(other.amse_forward),
145 amse_backward(other.amse_backward), N_data_forward(other.N_data_forward), N_data_backward(other.N_data_backward),
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),
147 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),
148 s(other.s), method(other.method), opt_method(other.opt_method), stepwise(other.stepwise), approximation(other.approximation),
149 num_models(other.num_models), seasonal(other.seasonal), stationary(other.stationary),
150 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) {
151 }
152
153 // Copy assignment operator
155 operator=(const InterpolARIMA &other) {
156 // protect against invalid self-assignment
157 if (this != &other) {
158 auto_arima_forward = auto_arima_copy(other.auto_arima_forward);
159 auto_arima_backward = auto_arima_copy(other.auto_arima_backward);
160
161 // 3: copy all the other fields from the other object
162 gap_loc = other.gap_loc;
163 N_gap = other.N_gap;
164 time = other.time;
165 pred_forward = other.pred_forward;
166 pred_backward = other.pred_backward;
167 norm = other.norm;
168 data = other.data;
169 xreg_vec_f = other.xreg_vec_f;
170 xreg_vec_b = other.xreg_vec_b;
171 data_forward = other.data_forward;
172 data_backward = other.data_backward;
173 new_xreg_vec_f = other.new_xreg_vec_f;
174 new_xreg_vec_b = other.new_xreg_vec_b;
175 N_data_forward = other.N_data_forward;
176 N_data_backward = other.N_data_backward;
177 max_p = other.max_p;
178 max_d = other.max_d;
179 max_q = other.max_q;
180 start_p = other.start_p;
181 start_q = other.start_q;
182 max_P = other.max_P;
183 max_D = other.max_D;
184 max_Q = other.max_Q;
185 start_P = other.start_P;
186 start_Q = other.start_Q;
187 r = other.r;
188 s = other.s;
189 method = other.method;
190 opt_method = other.opt_method;
191 stepwise = other.stepwise;
192 approximation = other.approximation;
193 num_models = other.num_models;
194 seasonal = other.seasonal;
195 stationary = other.stationary;
196
197 // 4: handle the pointers to the vectors
198 xreg_f = (xreg_vec_f.empty()) ? nullptr : &xreg_vec_f[0];
199 xreg_b = (xreg_vec_b.empty()) ? nullptr : &xreg_vec_b[0];
200 new_xreg_f = (new_xreg_vec_f.empty()) ? nullptr : &new_xreg_vec_f[0];
201 new_xreg_b = (new_xreg_vec_b.empty()) ? nullptr : &new_xreg_vec_b[0];
202 }
203 // by convention, always return *this
204 return *this;
205 }
206
208 auto_arima_free(auto_arima_forward);
209 auto_arima_free(auto_arima_backward);
210 delete xreg_f;
211 delete xreg_b;
212 delete new_xreg_f;
213 delete new_xreg_b;
214 }
215
216 // info
217 std::string toString();
218 std::string autoArimaInfo(auto_arima_object obj);
219
220
221private:
222 // Interpolation variables
223 Normalization norm;
224 std::vector<double> data;
225 size_t gap_loc;
226 size_t N_gap;
227 std::vector<double> time;
228 std::vector<double> pred_forward, pred_backward;
229
230 // Auto Arima variables
231 // const doesnt work with c
232 std::vector<double> xreg_vec_f, xreg_vec_b, data_forward, data_backward, new_xreg_vec_f, new_xreg_vec_b;
233 double *xreg_f;
234 double *xreg_b;
235 double *new_xreg_f;
236 double *new_xreg_b;
237 std::vector<double> amse_forward, amse_backward;
238 size_t N_data_forward, N_data_backward;
239 int max_p = 8, max_d = 3, max_q = 8;
240 int start_p = 2, start_q = 2;
241 int max_P = 2, max_D = 1, max_Q = 2;
242 int start_P = 1, start_Q = 1;
243 int r = 0, s = 0;
244 ObjectiveFunction method = CSS_MLE; OptimizationMethod opt_method = BFGS;
245 bool stepwise = true, approximation = true;
246 int num_models = 94;
247 bool seasonal = true, stationary = false;
248
249 bool consistencyCheck();
250 auto_arima_object initAutoArima(size_t N_data);
251
252 // last to be initialized
253public:
254 auto_arima_object auto_arima_forward;
255 auto_arima_object auto_arima_backward;
256
257private:
258 // (S)ARIMA variables
259 bool set_manual = false;
260 bool fill_backward_manual = false;
261
262public:
263 sarima_object sarima_forward;
264};
265
266} // namespace mio
267
268#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:263
InterpolARIMA & operator=(const InterpolARIMA &other)
Definition: InterpolARIMA.h:155
InterpolARIMA()
Definition: InterpolARIMA.cc:34
std::vector< double > getBackwardData()
Definition: InterpolARIMA.h:134
std::vector< double > getForwardData()
Definition: InterpolARIMA.h:133
void fillGap()
Definition: InterpolARIMA.cc:428
~InterpolARIMA()
Definition: InterpolARIMA.h:207
std::vector< double > simulate(int n_steps, int seed=0)
Definition: InterpolARIMA.cc:393
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:138
auto_arima_object auto_arima_forward
Definition: InterpolARIMA.h:254
std::vector< double > ARIMApredict(size_t n_steps)
Definition: InterpolARIMA.cc:401
std::vector< double > getData()
Definition: InterpolARIMA.h:132
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:543
void fillGapManual()
Definition: InterpolARIMA.cc:471
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:568
std::string autoArimaInfo(auto_arima_object obj)
Definition: InterpolARIMA.cc:175
void setVerbose(bool verbose=false)
Definition: InterpolARIMA.cc:370
auto_arima_object auto_arima_backward
Definition: InterpolARIMA.h:255
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