Loading [MathJax]/extensions/tex2jax.js
PISM, A Parallel Ice Sheet Model 2.2.1-cd005eec8 committed by Constantine Khrulev on 2025-03-07
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Interpolation1D.cc
Go to the documentation of this file.
1/* Copyright (C) 2015, 2016, 2017, 2018, 2019, 2021, 2023, 2024, 2025 PISM Authors
2 *
3 * This file is part of PISM.
4 *
5 * PISM is free software; you can redistribute it and/or modify it under the
6 * terms of the GNU General Public License as published by the Free Software
7 * Foundation; either version 3 of the License, or (at your option) any later
8 * version.
9 *
10 * PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12 * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13 * details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with PISM; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20#include <cstddef>
21#include <gsl/gsl_interp.h>
22#include <cassert>
23
24#include "pism/util/Interpolation1D.hh"
25#include "pism/util/error_handling.hh"
26
27namespace pism {
28
30 const std::vector<double> &input_x,
31 const std::vector<double> &output_x)
32 : Interpolation1D(type, input_x.data(), input_x.size(),
33 output_x.data(), output_x.size()) {
34 // empty
35}
36
38 const double *input_x, unsigned int input_x_size,
39 const double *output_x, unsigned int output_x_size) {
40
41 // the trivial case (the code below requires input_x_size >= 2)
42 if (input_x_size < 2) {
43 m_left.resize(output_x_size);
44 m_right.resize(output_x_size);
45 m_alpha.resize(output_x_size);
46
47 for (unsigned int k = 0; k < output_x_size; ++k) {
48 m_left[k] = 0.0;
49 m_right[k] = 0.0;
50 m_alpha[k] = 0.0;
51 }
52
53 return;
54 }
55
56 // input grid points have to be stored in the increasing order
57 for (unsigned int i = 0; i < input_x_size - 1; ++i) {
58 if (input_x[i] >= input_x[i + 1]) {
60 "an input grid for interpolation has to be "
61 "strictly increasing");
62 }
63 }
64
65 switch (type) {
66 case LINEAR:
67 init_linear(input_x, input_x_size, output_x, output_x_size);
68 break;
69 case NEAREST:
70 init_nearest(input_x, input_x_size, output_x, output_x_size);
71 break;
73 init_piecewise_constant(input_x, input_x_size, output_x, output_x_size);
74 break;
75 // LCOV_EXCL_START
76 default:
77 throw RuntimeError(PISM_ERROR_LOCATION, "invalid interpolation type");
78 // LCOV_EXCL_STOP
79 }
80}
81
82/**
83 * Compute linear interpolation indexes and weights.
84 *
85 * @param[in] input_x coordinates of the input grid
86 * @param[in] input_x_size number of points in the input grid
87 * @param[in] output_x coordinates of the output grid
88 * @param[in] output_x_size number of points in the output grid
89 */
90void Interpolation1D::init_linear(const double *input_x, unsigned int input_x_size,
91 const double *output_x, unsigned int output_x_size) {
92 assert(input_x_size >= 2);
93
94 m_left.resize(output_x_size);
95 m_right.resize(output_x_size);
96 m_alpha.resize(output_x_size);
97
98 // compute indexes and weights
99 for (unsigned int i = 0; i < output_x_size; ++i) {
100 double x = output_x[i];
101
102 // note: use "input_x_size" instead of "input_x_size - 1" to support extrapolation on
103 // the right
104 size_t
105 L = gsl_interp_bsearch(input_x, x, 0, input_x_size),
106 R = L + 1;
107
108 double alpha = 0.0;
109 if (x >= input_x[L] and R < input_x_size) {
110 // regular case
111 alpha = (x - input_x[L]) / (input_x[R] - input_x[L]);
112 } else {
113 // extrapolation
114 alpha = 0.0;
115 R = L;
116 }
117
118 assert(L < input_x_size);
119 assert(R < input_x_size);
120 assert(alpha >= 0.0 and alpha <= 1.0);
121
122 m_left[i] = static_cast<int>(L);
123 m_right[i] = static_cast<int>(R);
124 m_alpha[i] = alpha;
125 }
126}
127
128const std::vector<int>& Interpolation1D::left() const {
129 return m_left;
130}
131
132const std::vector<int>& Interpolation1D::right() const {
133 return m_right;
134}
135
136const std::vector<double>& Interpolation1D::alpha() const {
137 return m_alpha;
138}
139
140int Interpolation1D::left(size_t j) const {
141 return m_left[j];
142}
143
144int Interpolation1D::right(size_t j) const {
145 return m_right[j];
146}
147
148double Interpolation1D::alpha(size_t j) const {
149 return m_alpha[j];
150}
151
153 return (int)m_alpha.size();
154}
155
156std::vector<double> Interpolation1D::interpolate(const std::vector<double> &input_values) const {
157 std::vector<double> result(m_alpha.size());
158
159 interpolate(input_values.data(), result.data());
160
161 return result;
162}
163
164void Interpolation1D::interpolate(const double *input, double *output) const {
165 size_t n = m_alpha.size();
166 for (size_t k = 0; k < n; ++k) {
167 const int
168 L = m_left[k],
169 R = m_right[k];
170 output[k] = input[L] + m_alpha[k] * (input[R] - input[L]);
171 }
172}
173
174void Interpolation1D::init_nearest(const double *input_x, unsigned int input_x_size,
175 const double *output_x, unsigned int output_x_size) {
176
177 init_linear(input_x, input_x_size, output_x, output_x_size);
178
179 for (unsigned int j = 0; j < m_alpha.size(); ++j) {
180 m_alpha[j] = m_alpha[j] > 0.5 ? 1.0 : 0.0;
181 }
182}
183
184/*!
185 * Input grid `input_x` corresponds to *left* end-points of intervals.
186 */
187void Interpolation1D::init_piecewise_constant(const double *input_x, unsigned int input_x_size,
188 const double *output_x, unsigned int output_x_size) {
189 assert(input_x_size >= 2);
190
191 m_left.resize(output_x_size);
192 m_right.resize(output_x_size);
193 m_alpha.resize(output_x_size);
194
195 // compute indexes and weights
196 for (unsigned int i = 0; i < output_x_size; ++i) {
197
198 // note: use "input_x_size" instead of "input_x_size - 1" to support extrapolation on
199 // the right
200 size_t L = gsl_interp_bsearch(input_x, output_x[i], 0, input_x_size);
201
202 m_left[i] = static_cast<int>(L);
203 m_right[i] = static_cast<int>(L);
204 m_alpha[i] = 0.0;
205
206 assert(m_left[i] >= 0 and m_left[i] < (int)input_x_size);
207 assert(m_right[i] >= 0 and m_right[i] < (int)input_x_size);
208 assert(m_alpha[i] >= 0.0 and m_alpha[i] <= 1.0);
209 }
210}
211
212static std::map<size_t, double> weights_piecewise_constant(const double *x,
213 size_t x_size,
214 double a,
215 double b) {
216
217 size_t
218 al = gsl_interp_bsearch(x, a, 0, x_size),
219 ar = (a >= x[al] and al + 1 < x_size) ? al + 1 : al,
220 bl = gsl_interp_bsearch(x, b, 0, x_size),
221 br = (b >= x[bl] and bl + 1 < x_size) ? bl + 1 : bl;
222
223 std::map<size_t, double> result;
224
225 if (al == bl and ar == br) {
226 // both end points are in the same interval
227 result[al] += (b - a);
228 } else {
229 // first interval
230 result[al] += (x[ar] - a);
231
232 // intermediate intervals
233 for (size_t k = ar; k < bl; ++k) {
234 result[k] += x[k + 1] - x[k];
235 }
236
237 // last interval
238 result[bl] += b - x[bl];
239 }
240
241 return result;
242}
243
244static std::map<size_t, double> weights_piecewise_linear(const double *x,
245 size_t x_size,
246 double a,
247 double b) {
248
249 size_t
250 al = gsl_interp_bsearch(x, a, 0, x_size),
251 ar = (a >= x[al] and al + 1 < x_size) ? al + 1 : al,
252 bl = gsl_interp_bsearch(x, b, 0, x_size),
253 br = (b >= x[bl] and bl + 1 < x_size) ? bl + 1 : bl;
254
255 double
256 alpha_a = (al == ar) ? 0.0 : (a - x[al]) / (x[ar] - x[al]),
257 alpha_b = (bl == br) ? 0.0 : (b - x[bl]) / (x[br] - x[bl]);
258
259 std::map<size_t, double> result;
260
261 if (al == bl and ar == br) {
262 // both end points are in the same interval
263
264 result[al] += 0.5 * (2.0 - alpha_a - alpha_b) * (b - a);
265 result[ar] += 0.5 * (alpha_a + alpha_b) * (b - a);
266 } else {
267 // first interval
268 result[al] += 0.5 * (1.0 - alpha_a) * (x[ar] - a);
269 result[ar] += 0.5 * (1.0 + alpha_a) * (x[ar] - a);
270
271 // intermediate intervals
272 for (size_t k = ar; k < bl; ++k) {
273 size_t
274 L = k,
275 R = k + 1;
276 result[L] += 0.5 * (x[R] - x[L]);
277 result[R] += 0.5 * (x[R] - x[L]);
278 }
279
280 // last interval
281 result[bl] += 0.5 * (2.0 - alpha_b) * (b - x[bl]);
282 result[br] += 0.5 * alpha_b * (b - x[bl]);
283 }
284
285 return result;
286}
287
288/*!
289 * Compute weights for integrating a piece-wise linear or piece-wise constant function
290 * defined on the grid `x` from `a` to `b`.
291 *
292 * Uses constant extrapolation, both on the left and on the right.
293 *
294 * In the piece-wise constant case points in `x` are interpreted as left end points of
295 * intervals, i.e. the value `data[k]` corresponds to the interval `x[k], x[k + 1]`.
296 *
297 * To evaluate in the integral compute the dot product of data on the grid `x` with
298 * weights computed by this function:
299 *
300 * ```
301 * double result = 0.0;
302 * for (const auto &weight : weights) {
303 * size_t k = weight.first;
304 * double w = weight.second;
305 * result += w * data[k];
306 * }
307 * ```
308 */
309std::map<size_t, double> integration_weights(const double *x,
310 size_t x_size,
312 double a,
313 double b) {
314
315 if (a >= b) {
317 "invalid integration interval (a >= b)");
318 }
319
320 if (type != LINEAR and type != PIECEWISE_CONSTANT) {
322 "unsupported interpolation type");
323 }
324
325 auto weights = (type == LINEAR) ?
328
329 size_t N = x_size - 1;
330 double t0 = x[0];
331 double t1 = x[N];
332
333 // both points are to the left of [t0, t1]
334 if (b <= t0) {
335 return {{0, b - a}};
336 }
337
338 // both points are to the right of [t0, t1]
339 if (a >= t1) {
340 return {{N, b - a}};
341 }
342
343 // a is to the left of [t0, t1]
344 if (a < t0) {
345 auto W = weights(x, x_size, t0, b);
346 W[0] += t0 - a;
347 return W;
348 }
349
350 // b is to the right of [t0, t1]
351 if (b > t1) {
352 auto W = weights(x, x_size, a, t1);
353 W[N] += b - t1;
354 return W;
355 }
356
357 // both a and b are inside [t0, t1]
358 return weights(x, x_size, a, b);
359}
360
361std::map<size_t, double> integration_weights(const std::vector<double> &x,
363 double a,
364 double b) {
365 return integration_weights(x.data(), x.size(), type, a, b);
366}
367
368} // end of namespace pism
void init_linear(const double *input_x, unsigned int input_x_size, const double *output_x, unsigned int output_x_size)
void init_piecewise_constant(const double *input_x, unsigned int input_x_size, const double *output_x, unsigned int output_x_size)
Interpolation1D(InterpolationType type, const std::vector< double > &input_x, const std::vector< double > &output_x)
std::vector< int > m_right
std::vector< double > m_alpha
Interpolation weights.
const std::vector< double > & alpha() const
const std::vector< int > & right() const
void init_nearest(const double *input_x, unsigned int input_x_size, const double *output_x, unsigned int output_x_size)
const std::vector< int > & left() const
std::vector< int > m_left
Interpolation indexes.
std::vector< double > interpolate(const std::vector< double > &input_values) const
Return interpolated values (on the output grid) given input_values on the input grid.
#define PISM_ERROR_LOCATION
static const double L
Definition exactTestL.cc:40
#define n
Definition exactTestM.c:37
@ PIECEWISE_CONSTANT
std::map< size_t, double > integration_weights(const double *x, size_t x_size, InterpolationType type, double a, double b)
static std::map< size_t, double > weights_piecewise_linear(const double *x, size_t x_size, double a, double b)
static const double k
Definition exactTestP.cc:42
static std::map< size_t, double > weights_piecewise_constant(const double *x, size_t x_size, double a, double b)