MOONS Pipeline Reference Manual 0.13.2
moo_rebin.c
1/*
2 * This file is part of the MOONS Pipeline
3 * Copyright (C) 2002-2016 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19
20#ifdef HAVE_CONFIG_H
21#include <config.h>
22#endif
23
24/*-----------------------------------------------------------------------------
25 Includes
26 -----------------------------------------------------------------------------*/
27#include <math.h>
28#include <string.h>
29#include <cpl.h>
30#include <hdrl.h>
31#include "moo_params.h"
32#include "moo_single.h"
33#include "moo_badpix.h"
34#include "moo_ext.h"
35#include "moo_ext_single.h"
36#include "moo_utils.h"
37#include "moo_rebin.h"
38#include "moo_fits.h"
39#include "moo_pfits.h"
40#include "moo_qc.h"
41#include "moo_fibres_table.h"
42#ifdef _OPENMP
43#include <omp.h>
44#endif
45/*----------------------------------------------------------------------------*/
50/*----------------------------------------------------------------------------*/
53/*-----------------------------------------------------------------------------
54 Function codes
55 -----------------------------------------------------------------------------*/
56static cpl_table *
57_moo_rebin_fibre_cut_edges(cpl_image *data,
58 cpl_image *wdata,
59 cpl_image *errs,
60 cpl_image *qual,
61 int numfib,
62 int direction)
63{
64 cpl_table *result = NULL;
65 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
66 cpl_ensure(wdata != NULL, CPL_ERROR_NULL_INPUT, NULL);
67 cpl_ensure(errs != NULL, CPL_ERROR_NULL_INPUT, NULL);
68 cpl_ensure(qual != NULL, CPL_ERROR_NULL_INPUT, NULL);
69
70 cpl_vector *all_ref_flux = cpl_vector_new_from_image_row(data, numfib);
71 cpl_vector *all_ref_errs = cpl_vector_new_from_image_row(errs, numfib);
72 cpl_vector *all_ref_qual = cpl_vector_new_from_image_row(qual, numfib);
73 cpl_vector *all_ref_wave = cpl_vector_new_from_image_row(wdata, numfib);
74
76 int nx = cpl_image_get_size_x(data);
77 int ifirst = 0, ilast = nx - 1;
78 for (int i = 0; i < nx; i++) {
79 int q = (int)cpl_vector_get(all_ref_qual, i);
80 double w = cpl_vector_get(all_ref_wave, i);
81 if ((q & bp) == 0 && !isnan(w)) {
82 ifirst = i;
83 break;
84 }
85 }
86 for (int i = nx - 1; i >= 0; i--) {
87 int q = (int)cpl_vector_get(all_ref_qual, i);
88 double w = cpl_vector_get(all_ref_wave, i);
89 if ((q & bp) == 0 && !isnan(w)) {
90 ilast = i;
91 break;
92 }
93 }
94
95 int dsize = ilast - ifirst + 1;
96
97 result = cpl_table_new(dsize);
98 cpl_table_new_column(result, MOO_REBIN_REFTABLE_WAVE, CPL_TYPE_DOUBLE);
99 cpl_table_new_column(result, MOO_REBIN_REFTABLE_X, CPL_TYPE_DOUBLE);
100 cpl_table_new_column(result, MOO_REBIN_REFTABLE_Y, CPL_TYPE_INT);
101 cpl_table_new_column(result, MOO_REBIN_REFTABLE_FLUX, CPL_TYPE_DOUBLE);
102 cpl_table_new_column(result, MOO_REBIN_REFTABLE_ERR, CPL_TYPE_DOUBLE);
103 cpl_table_new_column(result, MOO_REBIN_REFTABLE_QUAL, CPL_TYPE_INT);
104
105 for (int i = ifirst; i <= ilast; i++) {
106 int it = i - ifirst;
107 double q = cpl_vector_get(all_ref_qual, i);
108 double f = cpl_vector_get(all_ref_flux, i);
109 double w = cpl_vector_get(all_ref_wave, i);
110 double e = cpl_vector_get(all_ref_errs, i);
111 int x = i + 1;
112 if (direction == -1) {
113 x = nx - i + 1;
114 }
115 cpl_table_set(result, MOO_REBIN_REFTABLE_WAVE, it, w);
116 cpl_table_set(result, MOO_REBIN_REFTABLE_X, it, x);
117 cpl_table_set(result, MOO_REBIN_REFTABLE_Y, it, numfib);
118 cpl_table_set(result, MOO_REBIN_REFTABLE_FLUX, it, f);
119 cpl_table_set(result, MOO_REBIN_REFTABLE_ERR, it, e);
120 cpl_table_set(result, MOO_REBIN_REFTABLE_QUAL, it, q);
121 }
122 cpl_propertylist *order = cpl_propertylist_new();
123 cpl_propertylist_append_bool(order, MOO_REBIN_REFTABLE_WAVE, CPL_FALSE);
124 cpl_table_sort(result, order);
125 cpl_propertylist_delete(order);
126 cpl_vector_delete(all_ref_wave);
127 cpl_vector_delete(all_ref_flux);
128 cpl_vector_delete(all_ref_errs);
129 cpl_vector_delete(all_ref_qual);
130
131 return result;
132}
133
134/* interpolate at the wavelength center of each bins */
135static cpl_error_code
136_moo_rebin_interpolate2(cpl_table *ref_table,
137 cpl_vector *res_pos,
138 cpl_vector *res_wave,
139 cpl_image *rbn_fluxes,
140 cpl_image *rbn_errs,
141 cpl_image *rbn_qual,
142 int imin,
143 int rbn_numfib)
144{
145 int ref_size = cpl_table_get_nrow(ref_table);
146 double *ref_wave_data =
147 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_WAVE);
148 double *ref_pos_data =
149 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_X);
150 int *ref_qual = cpl_table_get_data_int(ref_table, MOO_REBIN_REFTABLE_QUAL);
151 double *ref_flux_data =
152 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_FLUX);
153 double *ref_errs_data =
154 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_ERR);
155
156 cpl_vector *ref_wave = cpl_vector_wrap(ref_size, ref_wave_data);
157 cpl_vector *ref_pos = cpl_vector_wrap(ref_size, ref_pos_data);
158 cpl_vector *ref_flux = cpl_vector_wrap(ref_size, ref_flux_data);
159 cpl_vector *ref_errs = cpl_vector_wrap(ref_size, ref_errs_data);
160
161 /* second interpolation from xgrid to flux */
162 int res_size = cpl_vector_get_size(res_wave);
163 cpl_vector *res_flux = cpl_vector_new(res_size);
164 cpl_vector *res_errs = cpl_vector_new(res_size);
165 int *res_qual = cpl_calloc(res_size, sizeof(int));
166
167 cpl_bivector *ref_flux_points =
168 cpl_bivector_wrap_vectors(ref_pos, ref_flux);
169 cpl_bivector *res_points = cpl_bivector_wrap_vectors(res_pos, res_flux);
170 moo_interpolate_linear(res_points, res_errs, res_qual, ref_flux_points,
171 ref_errs, ref_qual);
172 cpl_bivector_unwrap_vectors(ref_flux_points);
173 cpl_bivector_unwrap_vectors(res_points);
174 cpl_vector_unwrap(ref_wave);
175 cpl_vector_unwrap(ref_pos);
176 cpl_vector_unwrap(ref_flux);
177 cpl_vector_unwrap(ref_errs);
178
179 for (int i = 0; i < res_size; i++) {
180 int rbn_idx = imin + 1 + i;
181 int q = res_qual[i];
182 double f = cpl_vector_get(res_flux, i);
183 double e = cpl_vector_get(res_errs, i);
184
185 cpl_image_set(rbn_fluxes, rbn_idx, rbn_numfib, f);
186 cpl_image_set(rbn_errs, rbn_idx, rbn_numfib, e);
187 cpl_image_set(rbn_qual, rbn_idx, rbn_numfib, q);
188 }
189#if MOO_DEBUG_REBIN
190 if (rbn_numfib == 379) {
191 cpl_table *t = cpl_table_new(res_size);
192 cpl_table_new_column(t, "WAVE", CPL_TYPE_DOUBLE);
193 cpl_table_new_column(t, "FLUX", CPL_TYPE_DOUBLE);
194 cpl_table_new_column(t, "FLUX_ERR", CPL_TYPE_DOUBLE);
195 cpl_table_new_column(t, "BPM", CPL_TYPE_INT);
196 for (int i = 0; i < res_size; i++) {
197 int rej;
198 double w = cpl_vector_get(res_wave, i);
199 double f = cpl_vector_get(res_flux, i);
200 double e = cpl_vector_get(res_errs, i);
201 int q = res_qual[i];
202 cpl_table_set_double(t, "WAVE", i, w);
203 cpl_table_set_double(t, "FLUX", i, f);
204 cpl_table_set_double(t, "FLUX_ERR", i, e);
205 cpl_table_set_int(t, "BPM", i, q);
206 }
207 cpl_table_save(t, NULL, NULL, "interpolate_linear_379_me.fits",
208 CPL_IO_CREATE);
209 cpl_table_delete(t);
210 }
211#endif
212 cpl_vector_delete(res_flux);
213 cpl_vector_delete(res_errs);
214 cpl_free(res_qual);
215
216 return CPL_ERROR_NONE;
217}
218
219static cpl_error_code
220_moo_rebin_sum(cpl_table *ref_table,
221 double step,
222 cpl_vector *res_pos,
223 cpl_vector *res_wave,
224 cpl_image *rbn_fluxes,
225 cpl_image *rbn_errs,
226 cpl_image *rbn_qual,
227 int imin,
228 int rbn_numfib)
229{
230 int ref_size = cpl_table_get_nrow(ref_table);
231 double *ref_wave_data =
232 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_WAVE);
233 double *ref_pos_data =
234 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_X);
235 int *ref_qual = cpl_table_get_data_int(ref_table, MOO_REBIN_REFTABLE_QUAL);
236 double *ref_flux_data =
237 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_FLUX);
238 double *ref_errs_data =
239 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_ERR);
240
241 cpl_vector *ref_wave = cpl_vector_wrap(ref_size, ref_wave_data);
242 cpl_vector *ref_start = cpl_vector_new(ref_size);
243 cpl_vector *ref_stop = cpl_vector_new(ref_size);
244 cpl_vector *ref_step = cpl_vector_new(ref_size);
245
246 double w1, w, w2, start, stop;
247 for (int i = 1; i < ref_size - 1; i++) {
248 w1 = cpl_vector_get(ref_wave, i - 1);
249 w = cpl_vector_get(ref_wave, i);
250 w2 = cpl_vector_get(ref_wave, i + 1);
251 start = w - (w - w1) / 2.;
252 stop = w + (w2 - w) / 2.;
253 cpl_vector_set(ref_start, i, start);
254 cpl_vector_set(ref_stop, i, stop);
255 cpl_vector_set(ref_step, i, step);
256 }
257
258 w = cpl_vector_get(ref_wave, 0);
259 w2 = cpl_vector_get(ref_wave, 1);
260 stop = w - (w2 - w) / 2.;
261 cpl_vector_set(ref_start, 0, w);
262 cpl_vector_set(ref_stop, 0, stop);
263 cpl_vector_set(ref_step, 0, step / 2.);
264
265 w1 = cpl_vector_get(ref_wave, ref_size - 2);
266 w = cpl_vector_get(ref_wave, ref_size - 1);
267 start = w - (w - w1) / 2.;
268 cpl_vector_set(ref_start, ref_size - 1, w1);
269 cpl_vector_set(ref_stop, ref_size - 1, w);
270 cpl_vector_set(ref_step, ref_size - 1, step / 2.);
271
272 cpl_vector *ref_pos = cpl_vector_wrap(ref_size, ref_pos_data);
273 cpl_vector *ref_flux = cpl_vector_wrap(ref_size, ref_flux_data);
274 cpl_vector *ref_errs = cpl_vector_wrap(ref_size, ref_errs_data);
275
276 int res_size = cpl_vector_get_size(res_pos);
277 int first = (int)cpl_vector_get(ref_pos, 0);
278
279
280 for (int i = 0; i < res_size - 1; i++) {
281 double d1 = cpl_vector_get(res_pos, i);
282 double d2 = cpl_vector_get(res_pos, i + 1);
283 double x1 = round(d1);
284 double x2 = round(d2);
285
286 int ix1 = (int)x1 - first;
287 int ix2 = (int)x2 - first;
288
289
290 double flux = 0.0;
291 double num = NAN;
292 double den = NAN;
293 double e = 0.0;
294 int q = MOO_BADPIX_GOOD;
295
296 if (ix1 == ix2) {
297 start = cpl_vector_get(res_wave, i);
298 stop = cpl_vector_get(res_wave, i + 1);
299 den = cpl_vector_get(ref_step, ix1);
300 num = (stop - start) / den;
301 e = cpl_vector_get(ref_errs, ix1) * num;
302 q |= ref_qual[ix1];
303 flux = cpl_vector_get(ref_flux, ix1) * num;
304 }
305 else {
306 start = cpl_vector_get(res_wave, i);
307 stop = cpl_vector_get(ref_stop, ix1);
308 den = cpl_vector_get(ref_step, ix1);
309 num = (stop - start) / den;
310 flux += cpl_vector_get(ref_flux, ix1) * num;
311 e = cpl_vector_get(ref_errs, ix1) * cpl_vector_get(ref_errs, ix1) *
312 num;
313 q |= ref_qual[ix1];
314
315 int min = ix1 + 1;
316 int max = ix2 - 1;
317 for (int j = min; j <= max; j++) {
318 q |= ref_qual[j];
319 start = cpl_vector_get(ref_start, j);
320 stop = cpl_vector_get(ref_stop, j);
321 den = cpl_vector_get(ref_step, j);
322 num = (stop - start) / den;
323 flux += cpl_vector_get(ref_flux, j) * num;
324 e += cpl_vector_get(ref_errs, j) * cpl_vector_get(ref_errs, j) *
325 num;
326 w = cpl_vector_get(ref_wave, j);
327 }
328 w = cpl_vector_get(ref_wave, ix2);
329 start = cpl_vector_get(ref_start, ix2);
330 stop = cpl_vector_get(res_wave, i + 1);
331 den = cpl_vector_get(ref_step, ix2);
332 num = (stop - start) / den;
333 flux += cpl_vector_get(ref_flux, ix2) * num;
334 e += cpl_vector_get(ref_errs, ix2) * cpl_vector_get(ref_errs, ix2) *
335 num;
336 q |= ref_qual[ix2];
337 e = sqrt(e);
338 }
339 int rbn_idx = imin + 1 + i;
340
341 cpl_image_set(rbn_fluxes, rbn_idx, rbn_numfib, flux);
342 cpl_image_set(rbn_errs, rbn_idx, rbn_numfib, e);
343 cpl_image_set(rbn_qual, rbn_idx, rbn_numfib, q);
344 }
345#if MOO_DEBUG_REBIN
346 if (rbn_numfib == 379) {
347 int nx = cpl_image_get_size_x(rbn_fluxes);
348 cpl_table *t = cpl_table_new(res_size - 1);
349 cpl_table_new_column(t, "WAVE", CPL_TYPE_DOUBLE);
350 cpl_table_new_column(t, "FLUX", CPL_TYPE_DOUBLE);
351 cpl_table_new_column(t, "FLUX_ERR", CPL_TYPE_DOUBLE);
352 cpl_table_new_column(t, "BPM", CPL_TYPE_INT);
353 for (int i = 0; i < res_size - 1; i++) {
354 int rej;
355 w1 = cpl_vector_get(res_wave, i);
356 w2 = cpl_vector_get(res_wave, i + 1);
357 w = (w1 + w2) / 2.;
358 int rbn_idx = imin + 1 + i;
359 double f = cpl_image_get(rbn_fluxes, rbn_idx, rbn_numfib, &rej);
360 double e = cpl_image_get(rbn_errs, rbn_idx, rbn_numfib, &rej);
361 int q = cpl_image_get(rbn_qual, rbn_idx, rbn_numfib, &rej);
362 cpl_table_set_double(t, "WAVE", i, w);
363 cpl_table_set_double(t, "FLUX", i, f);
364 cpl_table_set_double(t, "FLUX_ERR", i, e);
365 cpl_table_set_int(t, "BPM", i, q);
366 }
367 cpl_table_save(t, NULL, NULL, "integrate_linear_379_me.fits",
368 CPL_IO_CREATE);
369 cpl_table_delete(t);
370 }
371#endif
372 cpl_vector_unwrap(ref_pos);
373 cpl_vector_unwrap(ref_flux);
374 cpl_vector_unwrap(ref_errs);
375
376 cpl_vector_unwrap(ref_wave);
377 cpl_vector_delete(ref_start);
378 cpl_vector_delete(ref_stop);
379 cpl_vector_delete(ref_step);
380
381 return CPL_ERROR_NONE;
382}
383
384static cpl_error_code
385_moo_rebin_fibre(cpl_image *data,
386 cpl_image *wmap,
387 cpl_image *errs,
388 cpl_image *qual,
389 int numfib,
390 double step,
391 double min,
392 cpl_image *rbn_fluxes,
393 cpl_image *rbn_errs,
394 cpl_image *rbn_qual,
395 int fibfirst,
396 const char *method,
397 int conserv_flux,
398 int direction)
399{
400 cpl_table *ref_table = NULL;
401 cpl_vector *res_wave = NULL;
402 cpl_vector *res_pos = NULL;
403 cpl_bivector *ref_pos_points = NULL;
404 cpl_bivector *res_pos_points = NULL;
405 int rej;
406
407 int rbn_numfib = numfib + fibfirst;
408
409 ref_table =
410 _moo_rebin_fibre_cut_edges(data, wmap, errs, qual, numfib, direction);
411 int ref_size = cpl_table_get_nrow(ref_table);
412 double ref_wmin =
413 cpl_table_get_double(ref_table, MOO_REBIN_REFTABLE_WAVE, 0, &rej);
414 double ref_wmax = cpl_table_get_double(ref_table, MOO_REBIN_REFTABLE_WAVE,
415 ref_size - 1, &rej);
416
417 int rbn_size = cpl_image_get_size_x(rbn_fluxes);
418 double rbn_wmin = min;
419 double rbn_wmax = min + rbn_size * step;
420
421 if (ref_wmin < ref_wmax && ref_wmin < rbn_wmax) {
422 int imin = 0;
423 int imax = rbn_size;
424
425 imin = (int)ceil((ref_wmin - min) / step);
426 imax = floor((ref_wmax - min) / step);
427
428 if (imin < 0) {
429 imin = 0;
430 }
431 if (imin >= rbn_size) {
432 imin = rbn_size - 1;
433 }
434 if (imax < 0) {
435 imax = 0;
436 }
437 if (imax >= rbn_size) {
438 imax = rbn_size;
439 }
440
441 for (int i = 0; i < imin; i++) {
442 cpl_image_set(rbn_fluxes, i + 1, rbn_numfib, NAN);
443 cpl_image_set(rbn_errs, i + 1, rbn_numfib, NAN);
444 cpl_image_set(rbn_qual, i + 1, rbn_numfib,
446 }
447 for (int i = imax; i < rbn_size; i++) {
448 cpl_image_set(rbn_fluxes, i + 1, rbn_numfib, NAN);
449 cpl_image_set(rbn_errs, i + 1, rbn_numfib, NAN);
450 cpl_image_set(rbn_qual, i + 1, rbn_numfib,
452 }
453
454 if (strcmp(method, MOO_REBIN_METHOD_INTERPOLATE) == 0) {
455 /* first interpolate from lambda grid to xgrid */
456 int res_size = (imax - imin);
457 res_wave = cpl_vector_new(res_size);
458 res_pos = cpl_vector_new(res_size);
459
460 for (int i = 0; i < res_size; i++) {
461 double val1 = min + (imin + i) * step + step / 2.;
462 cpl_vector_set(res_wave, i, val1);
463 }
464 }
465 else {
466 /* first interpolate from lambda grid to xgrid */
467 int res_size = (imax - imin) + 1;
468 res_wave = cpl_vector_new(res_size);
469 res_pos = cpl_vector_new(res_size);
470 for (int i = 0; i < res_size; i++) {
471 double val1 = min + (imin + i) * step;
472 cpl_vector_set(res_wave, i, val1);
473 }
474 }
475
476 cpl_errorstate prev_state = cpl_errorstate_get();
477
478 double *ref_wave_data =
479 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_WAVE);
480 double *ref_pos_data =
481 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_X);
482 double *ref_flux_data =
483 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_FLUX);
484 double *ref_errs_data =
485 cpl_table_get_data_double(ref_table, MOO_REBIN_REFTABLE_ERR);
486
487 cpl_vector *ref_wave = cpl_vector_wrap(ref_size, ref_wave_data);
488 cpl_vector *ref_pos = cpl_vector_wrap(ref_size, ref_pos_data);
489 cpl_vector *ref_flux = cpl_vector_wrap(ref_size, ref_flux_data);
490 cpl_vector *ref_errs = cpl_vector_wrap(ref_size, ref_errs_data);
491 res_pos_points = cpl_bivector_wrap_vectors(res_wave, res_pos);
492 ref_pos_points = cpl_bivector_wrap_vectors(ref_wave, ref_pos);
493 cpl_bivector_interpolate_linear(res_pos_points, ref_pos_points);
494 cpl_bivector_unwrap_vectors(ref_pos_points);
495 cpl_bivector_unwrap_vectors(res_pos_points);
496 cpl_vector_unwrap(ref_wave);
497 cpl_vector_unwrap(ref_pos);
498 cpl_vector_unwrap(ref_flux);
499 cpl_vector_unwrap(ref_errs);
500
501 if (cpl_errorstate_is_equal(prev_state)) {
502 if (strcmp(method, MOO_REBIN_METHOD_INTERPOLATE) == 0) {
503 _moo_rebin_interpolate2(ref_table, res_pos, res_wave,
504 rbn_fluxes, rbn_errs, rbn_qual, imin,
505 rbn_numfib);
506 }
507 else {
508 _moo_rebin_sum(ref_table, step, res_pos, res_wave, rbn_fluxes,
509 rbn_errs, rbn_qual, imin, rbn_numfib);
510 }
511 }
512 else {
513 cpl_msg_error("moo_rebin", "#%d Invalid monotonicity (%f,%f)",
514 numfib, ref_wmin, ref_wmax);
515 cpl_errorstate_set(prev_state);
516 for (int i = 0; i < rbn_size; i++) {
517 cpl_image_set(rbn_fluxes, i + 1, rbn_numfib, NAN);
518 cpl_image_set(rbn_errs, i + 1, rbn_numfib, NAN);
519 cpl_image_set(rbn_qual, i + 1, rbn_numfib,
521 }
522 }
523
524 if (conserv_flux) {
525 /* first interpolate from lambda grid to xgrid */
526 int res_size = (imax - imin) + 1;
527 cpl_vector_delete(res_wave);
528 cpl_vector_delete(res_pos);
529 res_wave = cpl_vector_new(res_size);
530 res_pos = cpl_vector_new(res_size);
531 for (int i = 0; i < res_size; i++) {
532 double val1 = min + (imin + i) * step;
533 cpl_vector_set(res_wave, i, val1);
534 }
535 res_pos_points = cpl_bivector_wrap_vectors(res_wave, res_pos);
536 ref_wave = cpl_vector_wrap(ref_size, ref_wave_data);
537 ref_pos = cpl_vector_wrap(ref_size, ref_pos_data);
538 ref_pos_points = cpl_bivector_wrap_vectors(ref_wave, ref_pos);
539 cpl_bivector_interpolate_linear(res_pos_points, ref_pos_points);
540 for (int i = imin; i < imax; i++) {
541 int rbn_idx = i + 1;
542 double pix_lo = cpl_vector_get(res_pos, i - imin);
543 double pix_up = cpl_vector_get(res_pos, i - imin + 1);
544 double pix_size = pix_up - pix_lo;
545 double flux =
546 cpl_image_get(rbn_fluxes, rbn_idx, rbn_numfib, &rej);
547 cpl_image_set(rbn_fluxes, rbn_idx, rbn_numfib, flux * pix_size);
548 double err = cpl_image_get(rbn_errs, rbn_idx, rbn_numfib, &rej);
549 cpl_image_set(rbn_errs, rbn_idx, rbn_numfib, err * pix_size);
550 }
551 cpl_bivector_unwrap_vectors(ref_pos_points);
552 cpl_bivector_unwrap_vectors(res_pos_points);
553 cpl_vector_unwrap(ref_wave);
554 cpl_vector_unwrap(ref_pos);
555 }
556 }
557 else {
558 cpl_msg_error("moo_rebin",
559 "#%d Invalid solution (%f,%f) rbn range is '%f,%f)'",
560 numfib, ref_wmin, ref_wmax, rbn_wmin, rbn_wmax);
561 for (int i = 0; i < rbn_size; i++) {
562 cpl_image_set(rbn_fluxes, i + 1, rbn_numfib, NAN);
563 cpl_image_set(rbn_errs, i + 1, rbn_numfib, NAN);
564 cpl_image_set(rbn_qual, i + 1, rbn_numfib, MOO_BADPIX_CALIB_DEFECT);
565 }
566 }
567
568 cpl_vector_delete(res_wave);
569 cpl_vector_delete(res_pos);
570 cpl_table_delete(ref_table);
571 return CPL_ERROR_NONE;
572}
573
574/*----------------------------------------------------------------------------*/
593/*----------------------------------------------------------------------------*/
594static cpl_error_code
595_moo_rbn_compute_wmap_step(cpl_image *wmap,
596 double *qa,
597 double *qb,
598 int direction)
599{
600 int rej;
601 cpl_ensure_code(wmap != NULL, CPL_ERROR_NULL_INPUT);
602 int nx = cpl_image_get_size_x(wmap);
603 int ny = cpl_image_get_size_y(wmap);
604
605 cpl_image *pixsize = cpl_image_new(nx - 1, ny, CPL_TYPE_DOUBLE);
606 for (int j = 1; j <= ny; j++) {
607 for (int i = 1; i < nx - 1; i++) {
608 double diff = cpl_image_get(wmap, i + 1, j, &rej) -
609 cpl_image_get(wmap, i, j, &rej) * direction;
610 cpl_image_set(pixsize, i, j, diff);
611 }
612 }
613 moo_image_get_quartile(pixsize, qa, qb);
614
615 cpl_image_delete(pixsize);
616
617 return CPL_ERROR_NONE;
618}
619
620static cpl_error_code
621_moo_rebin_single(moo_rbn_single *rbn,
622 moo_ext_single *ext,
623 cpl_image *wmap,
624 const char *filename,
625 const int *health,
626 double step,
627 double min,
628 int first,
629 int *indexrbn,
630 const char **fibnames,
631 int *monotonous_tab,
632 cpl_array *indexes,
633 const char *method,
634 int conserv_flux,
635 int direction)
636{
637 double step_min, step_max;
638
639 cpl_ensure_code(rbn != NULL, CPL_ERROR_NULL_INPUT);
640 cpl_ensure_code(ext != NULL, CPL_ERROR_NULL_INPUT);
641 cpl_ensure_code(wmap != NULL, CPL_ERROR_NULL_INPUT);
642 cpl_ensure_code(filename != NULL, CPL_ERROR_NULL_INPUT);
643 cpl_ensure_code(health != NULL, CPL_ERROR_NULL_INPUT);
644
645 cpl_errorstate prestate = cpl_errorstate_get();
646 cpl_propertylist_copy_property_regexp(rbn->header, ext->header, "ESO DET *",
647 0);
648
649
650 _moo_rbn_compute_wmap_step(wmap, &step_min, &step_max, direction);
651
652 if ((strcmp(method, MOO_REBIN_METHOD_INTEGRATE) == 0) && step < step_min) {
653 cpl_msg_warning(
654 "moo_rebin",
655 "Use method %s with step %f but wavemap step is in range (%f,%f)",
656 method, step, step_min, step_max);
657 }
658 else if ((strcmp(method, MOO_REBIN_METHOD_INTERPOLATE) == 0) &&
659 step > step_max) {
660 cpl_msg_warning(
661 "moo_rebin",
662 "Use method %s with step %f but wavemap step is in range (%f,%f)",
663 method, step, step_min, step_max);
664 }
665
666 cpl_image *data = moo_ext_single_get_data(ext);
667 cpl_image *errs = moo_ext_single_get_errs(ext);
668 cpl_image *qual = moo_ext_single_get_qual(ext);
669
670 int nb_fibres = cpl_image_get_size_y(data);
671
672 cpl_image *rbn_fluxes = hdrl_image_get_image(rbn->image);
673 cpl_image *rbn_errs = hdrl_image_get_error(rbn->image);
674 cpl_image *rbn_qual = rbn->qual;
675#ifdef _OPENMP
676#pragma omp parallel default(none) \
677 shared(nb_fibres, health, rbn_fluxes, rbn_errs, rbn_qual, step, min, ext, \
678 data, wmap, errs, qual, first, indexes, method, indexrbn, \
679 fibnames, monotonous_tab, conserv_flux, direction)
680 {
681#pragma omp for
682#endif
683 for (int i = 1; i <= nb_fibres; i++) {
684 int idx = cpl_array_get_cplsize(indexes, i - 1, NULL);
685 int h = health[idx];
686 const char *fibname = fibnames[idx];
687
688 indexrbn[idx] = i + first;
689 if (h == 1) {
690 int monotonous = monotonous_tab[idx];
691 if (monotonous) {
692 _moo_rebin_fibre(data, wmap, errs, qual, i, step, min,
693 rbn_fluxes, rbn_errs, rbn_qual, first,
694 method, conserv_flux, direction);
695 }
696 else {
697 cpl_msg_warning(
698 __func__,
699 "skip #%s : no monotonous solution in WAVEMAP",
700 fibname);
701 }
702 }
703 }
704#ifdef _OPENMP
705 }
706#endif
707 // dump error from the state
708 if (!cpl_errorstate_is_equal(prestate)) {
709 cpl_msg_error(__func__, "Error in rebinning EXT:%s", ext->extname);
710 cpl_errorstate_dump(prestate, CPL_FALSE, cpl_errorstate_dump_one);
711 // Recover from the error(s) (Reset to prestate))
712 cpl_errorstate_set(prestate);
713 }
714 return CPL_ERROR_NONE;
715}
716
717static cpl_propertylist *
718_moo_create_primary_header(cpl_propertylist *ext_primary_header)
719{
720 cpl_propertylist *primary_header = NULL;
721
722 primary_header = cpl_propertylist_new();
723 cpl_propertylist_copy_property(primary_header, ext_primary_header,
724 MOO_PFITS_EXPTIME);
725 cpl_propertylist_copy_property(primary_header, ext_primary_header,
726 MOO_PFITS_MJDOBS);
727 cpl_propertylist_copy_property_regexp(primary_header, ext_primary_header,
728 "TM*", 0);
729 cpl_propertylist_copy_property_regexp(primary_header, ext_primary_header,
730 "ESO QC IS*", 0);
731 moo_pfits_set_fluxcal(primary_header, MOO_PFITS_FLUXCAL_UNCALIBRATED);
732 moo_qc_set_is_tellcor(primary_header, CPL_FALSE);
733 return primary_header;
734}
735
736static cpl_error_code
737_moo_update_single_header(cpl_propertylist *result_header,
738 int j,
739 double step,
740 double wlmin,
741 double wlmax)
742{
743 cpl_error_code status = CPL_ERROR_NONE;
744 moo_pfits_append_hduclass_data(result_header, j, -1);
745 cpl_propertylist_append_double(result_header, MOO_PFITS_CRPIX2,
746 MOO_EXT_SINGLE_CRPIX2);
747 cpl_propertylist_append_double(result_header, MOO_PFITS_CRVAL2,
748 MOO_EXT_SINGLE_CRVAL2);
749 cpl_propertylist_append_double(result_header, MOO_PFITS_CD1_2, 0.);
750 cpl_propertylist_append_double(result_header, MOO_PFITS_CD2_1, 0.);
751 cpl_propertylist_append_double(result_header, MOO_PFITS_CD2_2,
752 MOO_EXT_SINGLE_CDELT2);
753 cpl_propertylist_append_string(result_header, MOO_PFITS_CTYPE2,
754 MOO_EXT_SINGLE_CTYPE2);
755 cpl_propertylist_append_string(result_header, MOO_PFITS_BUNIT,
756 MOO_RBN_SINGLE_BUNIT);
757 moo_pfits_set_wlmin(result_header, wlmin);
758 moo_pfits_set_wlmax(result_header, wlmax);
759 moo_pfits_set_spec_bin(result_header, step);
760 return status;
761}
762/*----------------------------------------------------------------------------*/
781/*----------------------------------------------------------------------------*/
782moo_rbn *
783moo_rebin(moo_ext *ext,
784 moo_map *wmap,
785 moo_spectral_format *sformat,
786 moo_rebin_params *params,
787 const char *filename)
788{
789 moo_rbn *result = NULL;
790 const char *monotonous_colnames[] = { MOO_FIBRES_TABLE_MONOTONOUS_RI,
791 MOO_FIBRES_TABLE_MONOTONOUS_YJ,
792 MOO_FIBRES_TABLE_MONOTONOUS_H };
793
794 cpl_ensure(ext != NULL, CPL_ERROR_NULL_INPUT, NULL);
795 cpl_ensure(wmap != NULL, CPL_ERROR_NULL_INPUT, NULL);
796 cpl_ensure(sformat != NULL, CPL_ERROR_NULL_INPUT, NULL);
797 cpl_ensure(params != NULL, CPL_ERROR_NULL_INPUT, NULL);
798 const char *method = params->method;
799 cpl_ensure(method != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
800
801 cpl_errorstate prestate = cpl_errorstate_get();
802 unsigned int badpix_level = MOO_BADPIX_GOOD;
803
804 cpl_table *fibres_table = moo_ext_get_fibre_table(ext);
805 cpl_ensure(fibres_table != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
806
807 cpl_msg_info(__func__, "Rebin spectra using method %s flux conservation:%d",
808 method, params->conserv_flux);
809
810
811 moo_fits_create(filename);
812
813 result = moo_rbn_new();
814
815 result->filename = filename;
816
817 result->primary_header = _moo_create_primary_header(ext->primary_header);
818 cpl_image **maptab = wmap->data;
819 cpl_propertylist **headertab = wmap->data_header;
820 cpl_table *wmap_fibre_table = wmap->fibre_table;
821 cpl_table *rbn_fibre_table = cpl_table_duplicate(fibres_table);
822 moo_fibres_table_add_rebin_cols(rbn_fibre_table);
823
824 result->fibre_table = rbn_fibre_table;
825 cpl_ensure(rbn_fibre_table != NULL, CPL_ERROR_NULL_INPUT, NULL);
826 cpl_array *idx1 = NULL, *idx2 = NULL;
827
828 moo_try_check(
829 (idx1 = moo_fibres_table_get_spectro_indexext(rbn_fibre_table, 1),
830 idx2 = moo_fibres_table_get_spectro_indexext(rbn_fibre_table, 2)),
831 " Can't get index from INDEXEXT");
832
833
834 const int *health =
835 cpl_table_get_data_int_const(rbn_fibre_table, MOO_FIBRES_TABLE_HEALTH);
836 int *indexrbn =
837 cpl_table_get_data_int(rbn_fibre_table, MOO_FIBRES_TABLE_INDEXRBN);
838 const char **fibnames =
839 cpl_table_get_data_string_const(rbn_fibre_table,
840 MOO_FIBRES_TABLE_FIBRE);
841
842 int nrows = cpl_table_get_nrow(rbn_fibre_table);
843 for (int i = 0; i < nrows; i++) {
844 const char *fibname = fibnames[i];
845 cpl_table_select_all(wmap_fibre_table);
847 MOO_FIBRES_TABLE_FIBRE, fibname);
848 cpl_array *sel = cpl_table_where_selected(wmap_fibre_table);
849
850 int idx = cpl_array_get_cplsize(sel, 0, NULL);
851 for (int j = 0; j < 3; j++) {
852 if (cpl_table_has_column(wmap_fibre_table,
853 monotonous_colnames[j])) {
854 int res = cpl_table_get_int(wmap_fibre_table,
855 monotonous_colnames[j], idx, NULL);
856 cpl_table_set_int(rbn_fibre_table, monotonous_colnames[j], i,
857 res);
858 }
859 else {
860 cpl_table_set_int(rbn_fibre_table, monotonous_colnames[j], i,
861 1);
862 }
863 }
864 cpl_array_delete(sel);
865 }
866
867 cpl_array *idx_tab[] = { idx1, idx2 };
868
869 int size1 = cpl_array_get_size(idx1);
870 int first_tab[] = { 0, size1 };
871
872 cpl_msg_indent_more();
873
874 for (int j = 0; j < 3; j++) {
875 double wmin, wmax;
876 moo_spectral_format_get_wave_range(sformat, j, &wmin, &wmax);
877 double step = params->step[j];
878 moo_rbn_single *rbn_single = moo_rbn_single_new(j);
879 int nx = floor((wmax - wmin) / step);
880 int ny = cpl_table_get_nrow(rbn_fibre_table);
881 int *monotonous_tab =
882 cpl_table_get_data_int(rbn_fibre_table, monotonous_colnames[j]);
883 _moo_update_single_header(rbn_single->header, j, step, wmin, wmax);
884 moo_rbn_single_set_wcs1(rbn_single, 0.5, wmin, step, MOO_REBIN_CTYPE1,
885 MOO_REBIN_CUNIT1);
886 for (int i = 1; i <= 2; i++) {
887 cpl_array *idx = idx_tab[i - 1];
888 int first = first_tab[i - 1];
889 cpl_image *map = maptab[(i - 1) * 3 + j];
890 cpl_propertylist *map_header = headertab[(i - 1) * 3 + j];
891 moo_ext_single *ext_single =
892 moo_ext_load_single(ext, j, i, badpix_level);
893 if (ext_single != NULL && map != NULL) {
894 if (i == 1) {
895 cpl_msg_info(
896 "moo_rebin",
897 "Rebin extension %s with step %f and range [%f,%f]",
898 moo_detector_get_name(j), step, wmin, wmax);
899 }
900 if (rbn_single->image == NULL) {
901 rbn_single->image = hdrl_image_new(nx, ny);
902 rbn_single->qual = cpl_image_new(nx, ny, CPL_TYPE_INT);
903
904 cpl_image *img = moo_rbn_single_get_data(rbn_single);
905 cpl_image *error = moo_rbn_single_get_errs(rbn_single);
906 double *idata = cpl_image_get_data(img);
907 double *edata = cpl_image_get_data(error);
908 int *qdata = cpl_image_get_data(rbn_single->qual);
909 for (int k = 0; k < nx * ny; k++) {
911 edata[k] = NAN;
912 idata[k] = NAN;
913 }
914 }
915 cpl_propertylist_copy_property(rbn_single->header,
916 ext_single->header,
917 MOO_PFITS_BUNIT);
918 moo_spectral_format_info *info =
919 moo_spectral_format_get(sformat, j, i);
920 if (cpl_propertylist_has(
921 map_header, MOONS_QC_WAVECAL_MONOTONOUS_SOLUTION)) {
922 int monotonous = cpl_propertylist_get_bool(
923 map_header, MOONS_QC_WAVECAL_MONOTONOUS_SOLUTION);
924 if (!monotonous) {
925 cpl_msg_warning(__func__,
926 "WAVEMAP %s: not monotonous solution "
927 "for all fibres",
928 ext_single->extname);
929 }
930 }
931
932 _moo_rebin_single(rbn_single, ext_single, map, filename, health,
933 step, wmin, first, indexrbn, fibnames,
934 monotonous_tab, idx, method,
935 params->conserv_flux, info->direction);
936
937
939 }
940 }
941 if (rbn_single->image == NULL) {
942 moo_rbn_single_delete(rbn_single);
943 rbn_single = NULL;
944 }
945 moo_rbn_add_single(result, j, rbn_single);
946 }
947 moo_rbn_add_fibre_table(result, rbn_fibre_table);
948
949 cpl_msg_indent_less();
950moo_try_cleanup:
951 cpl_array_delete(idx1);
952 cpl_array_delete(idx2);
953
954 if (!cpl_errorstate_is_equal(prestate)) {
955 moo_rbn_delete(result);
956 result = NULL;
957 cpl_msg_error(__func__, "Error in rebin for EXT %s", ext->filename);
958 cpl_errorstate_dump(prestate, CPL_FALSE, cpl_errorstate_dump_one);
959 // Recover from the error(s) (Reset to prestate))
960 cpl_errorstate_set(prestate);
961 }
962 return result;
963}
#define MOO_BADPIX_OUTSIDE_DATA_RANGE
Definition: moo_badpix.h:62
#define MOO_BADPIX_CALIB_DEFECT
Definition: moo_badpix.h:48
#define MOO_BADPIX_GOOD
Definition: moo_badpix.h:36
const char * moo_detector_get_name(moo_detector_type type)
Get the extension name of a detector.
Definition: moo_detector.c:183
cpl_image * moo_ext_single_get_qual(moo_ext_single *self)
Get image of qual.
cpl_image * moo_ext_single_get_errs(moo_ext_single *self)
Get image of errs.
cpl_image * moo_ext_single_get_data(moo_ext_single *self)
Get image of data.
moo_ext_single * moo_ext_load_single(moo_ext *self, moo_detector_type type, int num, unsigned int level)
Load the type part in EXT and return it.
Definition: moo_ext.c:156
cpl_table * moo_ext_get_fibre_table(moo_ext *self)
Get the FIBRE TABLE in EXT.
Definition: moo_ext.c:344
cpl_array * moo_fibres_table_get_spectro_indexext(cpl_table *table, int num)
get the index of a spectro in the fibre table sort by spetcro,indexext
cpl_error_code moo_fibres_table_add_rebin_cols(cpl_table *table)
add rebin additional columns
cpl_error_code moo_fits_create(const char *filename)
Create a new fits file with empty propertylist.
Definition: moo_fits.c:533
moo_rbn_single * moo_rbn_single_new(moo_detector_type type)
Create a new moo_rbn_single.
void moo_rbn_single_delete(moo_rbn_single *self)
Delete a moo_rbn_single.
cpl_image * moo_rbn_single_get_data(moo_rbn_single *self)
Get image of data.
cpl_error_code moo_rbn_single_set_wcs1(moo_rbn_single *self, double crpix1, double crval1, double cd1_1, const char *ctype1, const char *cunit1)
Set the WCS1 of the extension.
cpl_image * moo_rbn_single_get_errs(moo_rbn_single *self)
Get image of errs.
void moo_rbn_delete(moo_rbn *self)
Delete a moo_rbn.
Definition: moo_rbn.c:120
cpl_error_code moo_rbn_add_single(moo_rbn *self, moo_detector_type type, moo_rbn_single *single)
Add RBN_SINGLE extension to RBN filename and update moo_rbn structure.
Definition: moo_rbn.c:226
cpl_error_code moo_rbn_add_fibre_table(moo_rbn *self, cpl_table *fibre_table)
Add fibre table to RBN filename and update moo_rbn structure.
Definition: moo_rbn.c:288
moo_rbn * moo_rbn_new(void)
Create a new moo_rbn.
Definition: moo_rbn.c:67
void moo_spectral_format_info_delete(moo_spectral_format_info *self)
Delete a moo_spectral_format_info.
moo_spectral_format_info * moo_spectral_format_get(moo_spectral_format *self, moo_detector_type type, int ntas)
Get the spectral format information for a given ARM.
cpl_error_code moo_spectral_format_get_wave_range(moo_spectral_format *self, moo_detector_type type, double *min, double *max)
Get the spectral format wave range for a given detector.
moo_rbn * moo_rebin(moo_ext *ext, moo_map *wmap, moo_spectral_format *sformat, moo_rebin_params *params, const char *filename)
Rebin an EXT spectra in RBN format.
Definition: moo_rebin.c:783
cpl_error_code moo_pfits_set_spec_bin(cpl_propertylist *plist, double val)
Set the SPEC_BIN Keyword.
Definition: moo_pfits.c:1718
cpl_error_code moo_pfits_set_wlmin(cpl_propertylist *plist, double val)
Set the WAVELMIN Keyword.
Definition: moo_pfits.c:1666
cpl_error_code moo_pfits_append_hduclass_data(cpl_propertylist *plist, moo_detector_type type, int ntas)
Set the HDUCLASS DATA Keyword.
Definition: moo_pfits.c:1504
cpl_error_code moo_pfits_set_wlmax(cpl_propertylist *plist, double val)
Set the WAVELMAX Keyword.
Definition: moo_pfits.c:1692
cpl_error_code moo_pfits_set_fluxcal(cpl_propertylist *plist, const char *val)
Set the FLUXCAL Keyword.
Definition: moo_pfits.c:1640
cpl_error_code moo_qc_set_is_tellcor(cpl_propertylist *plist, cpl_boolean val)
Set the QC.IS.TELLCOR value.
Definition: moo_qc.c:2451
cpl_size moo_table_and_selected_sequal_string(cpl_table *table, const char *name, const char *string)
Select from unselected table rows, by comparing column values with a constant.
Definition: moo_utils.c:2317
cpl_error_code moo_interpolate_linear(cpl_bivector *fout, cpl_vector *fout_errs, int *fout_qual, const cpl_bivector *fref, const cpl_vector *fref_errs, const int *fref_qual)
Linear interpolation of a 1d-function.
Definition: moo_utils.c:1160
cpl_error_code moo_image_get_quartile(cpl_image *image, double *qmin, double *qmax)
compute first and last quartile from an image
Definition: moo_utils.c:1906